Project
Loading...
Searching...
No Matches
o2::math_utils Namespace Reference

Namespaces

namespace  detail
 

Classes

class  CachingTF1
 
class  Chebyshev3D
 
class  Chebyshev3DCalc
 
class  Legendre1DPolynominal
 
class  Legendre2DPolynominal
 .... More...
 
class  RandomRing
 
class  Rotation2D
 transformation types More...
 
struct  StatisticsData
 
class  SymMatrixSolver
 
class  Transform3D
 
struct  TransformType
 
struct  Tsallis
 

Typedefs

using Rotation2Df_t = Rotation2D< float >
 
using Rotation2Dd_t = Rotation2D< double >
 
template<typename T >
using CircleXY = detail::CircleXY< T >
 
using CircleXYf_t = detail::CircleXY< float >
 
using CircleXYd_t = detail::CircleXY< double >
 
template<typename T >
using IntervalXY = detail::IntervalXY< T >
 
using IntervalXYf_t = detail::IntervalXY< float >
 
using IntervalXYd_t = detail::IntervalXY< double >
 
template<typename T >
using Bracket = detail::Bracket< T >
 
using Bracketf_t = detail::Bracket< float >
 
using Bracketd_t = detail::Bracket< double >
 
template<typename T >
using Point2D = ROOT::Math::PositionVector2D< ROOT::Math::Cartesian2D< T >, ROOT::Math::DefaultCoordinateSystemTag >
 
template<typename T >
using Vector2D = ROOT::Math::DisplacementVector2D< ROOT::Math::Cartesian2D< T >, ROOT::Math::DefaultCoordinateSystemTag >
 
template<typename T >
using Point3D = ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< T >, ROOT::Math::DefaultCoordinateSystemTag >
 
template<typename T >
using Vector3D = ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< T >, ROOT::Math::DefaultCoordinateSystemTag >
 
template<typename T , uint32_t N>
using SVector = ROOT::Math::SVector< T, N >
 
template<class T , uint32_t D1, uint32_t D2, class R >
using SMatrix = ROOT::Math::SMatrix< T, D1, D2, R >
 
template<class T , uint32_t D>
using MatRepSym = ROOT::Math::MatRepSym< T, D >
 
template<class T , uint32_t D1, uint32_t D2 = D1>
using MatRepStd = ROOT::Math::MatRepStd< T, D1, D2 >
 

Functions

template<class T , unsigned int D>
Dot (const SVector< T, D > &lhs, const SVector< T, D > &rhs)
 
template<class T , unsigned int D1, unsigned int D2, class R >
SMatrix< T, D1, D1, MatRepSym< T, D1 > > Similarity (const SMatrix< T, D1, D2, R > &lhs, const SMatrix< T, D2, D2, MatRepSym< T, D2 > > &rhs)
 
template<typename T >
TFitResultPtr fit (const size_t nBins, const T *arr, const T xMin, const T xMax, TF1 &func, std::string_view option="")
 
template<typename T >
bool medmadGaus (size_t nBins, const T *arr, const T xMin, const T xMax, std::array< double, 3 > &param)
 
template<typename T >
Double_t fitGaus (const size_t nBins, const T *arr, const T xMin, const T xMax, std::vector< T > &param)
 
template<typename T >
double fitGaus (size_t nBins, const T *arr, const T xMin, const T xMax, std::array< double, 3 > &param, ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > *covMat=nullptr, int minVal=2, bool applyMAD=true)
 
template<typename T >
StatisticsData getStatisticsData (const T *arr, const size_t nBins, const double xMin, const double xMax)
 
template<typename T , typename R = double>
R median (std::vector< T > v)
 
template<typename T >
void SortData (std::vector< T > const &values, std::vector< size_t > &index)
 
template<typename T >
bool LTMUnbinned (const std::vector< T > &data, std::vector< size_t > &index, std::array< float, 7 > &params, float fracKeep)
 
template<typename T >
void Reorder (std::vector< T > &data, const std::vector< size_t > &index)
 
template<typename T >
bool LTMUnbinnedSig (const std::vector< T > &data, std::vector< size_t > &index, std::array< float, 7 > &params, float fracKeepMin, float sigTgt, bool sorted=false)
 
template<typename T >
selKthMin (int k, int np, T *arr)
 
template<typename T >
MAD2Sigma (int np, T *y)
 
 GPUdi () float to02Pi(float phi)
 
 GPUdi () void bringTo02Pi(float &phi)
 
float toPMPiGen (float phi)
 
double toPMPiGend (double phi)
 
void bringToPMPiGen (float &phi)
 
void bringToPMPiGend (double &phi)
 
float to02PiGen (float phi)
 
double to02PiGend (double phi)
 
void bringTo02PiGen (float &phi)
 
void bringTo02PiGend (double &phi)
 
float toPMPi (float phi)
 
double toPMPid (double phi)
 
void bringToPMPi (float &phi)
 
void bringToPMPid (double &phi)
 
std::tuple< float, float > rotateZInv (float xG, float yG, float snAlp, float csAlp)
 
std::tuple< double, double > rotateZInvd (double xG, double yG, double snAlp, double csAlp)
 
std::tuple< float, float > rotateZ (float xL, float yL, float snAlp, float csAlp)
 
std::tuple< double, double > rotateZd (double xL, double yL, double snAlp, double csAlp)
 
void rotateZ (std::array< float, 3 > &xy, float alpha)
 
void rotateZd (std::array< double, 3 > &xy, double alpha)
 
int angle2Sector (float phi)
 
int angle2Sectord (double phi)
 
float sector2Angle (int sect)
 
double sector2Angled (int sect)
 
float angle2Alpha (float phi)
 
double angle2Alphad (double phi)
 
 GPUhdi () float fastATan2(float y
 
template<class T >
 GPUhdi () T min(const T x
 

Variables

float & s
 
float float & c
 
float yL
 
float float & xG
 
float float float & yG
 
float float float float snAlp
 
float float float float float csAlp
 
float float & xL
 
float x
 
const T y
 

Typedef Documentation

◆ Bracket

template<typename T >
using o2::math_utils::Bracket = typedef detail::Bracket<T>

Definition at line 39 of file Primitive2D.h.

◆ Bracketd_t

Definition at line 41 of file Primitive2D.h.

◆ Bracketf_t

Definition at line 40 of file Primitive2D.h.

◆ CircleXY

template<typename T >
using o2::math_utils::CircleXY = typedef detail::CircleXY<T>

Definition at line 29 of file Primitive2D.h.

◆ CircleXYd_t

Definition at line 31 of file Primitive2D.h.

◆ CircleXYf_t

Definition at line 30 of file Primitive2D.h.

◆ IntervalXY

template<typename T >
using o2::math_utils::IntervalXY = typedef detail::IntervalXY<T>

Definition at line 34 of file Primitive2D.h.

◆ IntervalXYd_t

Definition at line 36 of file Primitive2D.h.

◆ IntervalXYf_t

Definition at line 35 of file Primitive2D.h.

◆ MatRepStd

template<class T , uint32_t D1, uint32_t D2 = D1>
using o2::math_utils::MatRepStd = typedef ROOT::Math::MatRepStd<T, D1, D2>

Definition at line 63 of file GPUROOTSMatrixFwd.h.

◆ MatRepSym

template<class T , uint32_t D>
using o2::math_utils::MatRepSym = typedef ROOT::Math::MatRepSym<T, D>

Definition at line 61 of file GPUROOTSMatrixFwd.h.

◆ Point2D

template<typename T >
using o2::math_utils::Point2D = typedef ROOT::Math::PositionVector2D<ROOT::Math::Cartesian2D<T>, ROOT::Math::DefaultCoordinateSystemTag>

Definition at line 64 of file GPUROOTCartesianFwd.h.

◆ Point3D

template<typename T >
using o2::math_utils::Point3D = typedef ROOT::Math::PositionVector3D<ROOT::Math::Cartesian3D<T>, ROOT::Math::DefaultCoordinateSystemTag>

Definition at line 68 of file GPUROOTCartesianFwd.h.

◆ Rotation2Dd_t

using o2::math_utils::Rotation2Dd_t = typedef Rotation2D<double>

Definition at line 155 of file Cartesian.h.

◆ Rotation2Df_t

using o2::math_utils::Rotation2Df_t = typedef Rotation2D<float>

Definition at line 154 of file Cartesian.h.

◆ SMatrix

template<class T , uint32_t D1, uint32_t D2, class R >
using o2::math_utils::SMatrix = typedef ROOT::Math::SMatrix<T, D1, D2, R>

Definition at line 59 of file GPUROOTSMatrixFwd.h.

◆ SVector

template<typename T , uint32_t N>
using o2::math_utils::SVector = typedef ROOT::Math::SVector<T, N>

Definition at line 57 of file GPUROOTSMatrixFwd.h.

◆ Vector2D

template<typename T >
using o2::math_utils::Vector2D = typedef ROOT::Math::DisplacementVector2D<ROOT::Math::Cartesian2D<T>, ROOT::Math::DefaultCoordinateSystemTag>

Definition at line 66 of file GPUROOTCartesianFwd.h.

◆ Vector3D

template<typename T >
using o2::math_utils::Vector3D = typedef ROOT::Math::DisplacementVector3D<ROOT::Math::Cartesian3D<T>, ROOT::Math::DefaultCoordinateSystemTag>

Definition at line 70 of file GPUROOTCartesianFwd.h.

Function Documentation

◆ angle2Alpha()

float o2::math_utils::angle2Alpha ( float  phi)
inline

Definition at line 203 of file Utils.h.

◆ angle2Alphad()

double o2::math_utils::angle2Alphad ( double  phi)
inline

Definition at line 208 of file Utils.h.

◆ angle2Sector()

int o2::math_utils::angle2Sector ( float  phi)
inline

Definition at line 183 of file Utils.h.

◆ angle2Sectord()

int o2::math_utils::angle2Sectord ( double  phi)
inline

Definition at line 188 of file Utils.h.

◆ bringTo02PiGen()

void o2::math_utils::bringTo02PiGen ( float &  phi)
inline

Definition at line 80 of file Utils.h.

◆ bringTo02PiGend()

void o2::math_utils::bringTo02PiGend ( double &  phi)
inline

Definition at line 85 of file Utils.h.

◆ bringToPMPi()

void o2::math_utils::bringToPMPi ( float &  phi)
inline

Definition at line 100 of file Utils.h.

◆ bringToPMPid()

void o2::math_utils::bringToPMPid ( double &  phi)
inline

Definition at line 105 of file Utils.h.

◆ bringToPMPiGen()

void o2::math_utils::bringToPMPiGen ( float &  phi)
inline

Definition at line 60 of file Utils.h.

◆ bringToPMPiGend()

void o2::math_utils::bringToPMPiGend ( double &  phi)
inline

Definition at line 65 of file Utils.h.

◆ Dot()

template<class T , unsigned int D>
T o2::math_utils::Dot ( const SVector< T, D > &  lhs,
const SVector< T, D > &  rhs 
)
inline

Definition at line 257 of file Cartesian.h.

◆ fit()

template<typename T >
TFitResultPtr o2::math_utils::fit ( const size_t  nBins,
const T *  arr,
const T  xMin,
const T  xMax,
TF1 func,
std::string_view  option = "" 
)

fit 1D array of histogrammed data with generic root function

The code was extracted out of ROOT to be able to do fitting on an array with histogrammed data instead of root histograms. It is a stripped down version, so does not provide the same functionality. To be used with care.

Parameters
[in]nbinssize of the array and number of histogram bins
[in]arrarray with elements
[in]xMinminimum range of the array
[in]xMaxmaximum range of the array
[in]funcfit function

Definition at line 59 of file fit.h.

◆ fitGaus() [1/2]

template<typename T >
Double_t o2::math_utils::fitGaus ( const size_t  nBins,
const T *  arr,
const T  xMin,
const T  xMax,
std::vector< T > &  param 
)

fast fit of an array with ranges (histogram) with gaussian function

Fitting procedure:

  1. Step - make logarithm
  2. Linear fit (parabola) - more robust, always converges, fast
Parameters
[in]nbinssize of the array and number of histogram bins
[in]arrarray with elements
[in]xMinminimum range of the array
[in]xMaxmaximum range of the array
[out]paramreturn paramters of the fit (0-Constant, 1-Mean, 2-Sigma, 3-Sum)
Returns
chi2 or exit code >0: the chi2 returned by TLinearFitter -3: only three points have been used for the calculation - no fitter was used -2: only two points have been used for the calculation - center of gravity was uesed for calculation -1: only one point has been used for the calculation - center of gravity was uesed for calculation -4: invalid result!!

Definition at line 231 of file fit.h.

◆ fitGaus() [2/2]

template<typename T >
double o2::math_utils::fitGaus ( size_t  nBins,
const T *  arr,
const T  xMin,
const T  xMax,
std::array< double, 3 > &  param,
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > *  covMat = nullptr,
int  minVal = 2,
bool  applyMAD = true 
)

Definition at line 359 of file fit.h.

◆ getStatisticsData()

template<typename T >
StatisticsData o2::math_utils::getStatisticsData ( const T *  arr,
const size_t  nBins,
const double  xMin,
const double  xMax 
)

calculate statistical parameters on a binned array

The function assumes a binned array of

Parameters
nBinssize of the array
xMinlower histogram bound
xMaxupper histogram bound
Todo:
make return type templated?

Definition at line 470 of file fit.h.

◆ GPUdi() [1/2]

o2::math_utils::GPUdi ( )

Definition at line 30 of file Utils.h.

◆ GPUdi() [2/2]

o2::math_utils::GPUdi ( ) &

Definition at line 40 of file Utils.h.

◆ GPUhdi() [1/2]

o2::math_utils::GPUhdi ( )

Definition at line 245 of file Utils.h.

◆ GPUhdi() [2/2]

template<class T >
o2::math_utils::GPUhdi ( ) const

◆ LTMUnbinned()

template<typename T >
bool o2::math_utils::LTMUnbinned ( const std::vector< T > &  data,
std::vector< size_t > &  index,
std::array< float, 7 > &  params,
float  fracKeep 
)

LTM : Trimmed mean of unbinned array

Robust statistic to estimate properties of the distribution To handle binning error special treatment for definition of unbinned data see: http://en.wikipedia.org/w/index.php?title=Trimmed_estimator&oldid=582847999

Parameters
dataInput vector (unsorted)
indexVector with indices of sorted input data
paramsArray with characteristics of distribution
fracKeepFraction of data to be kept
Returns
Flag if successfull Characteristics:
  1. area
  2. mean
  3. rms
  4. error estimate of mean
  5. error estimate of RMS
  6. first accepted element (of sorted array)
  7. last accepted element (of sorted array)

Definition at line 569 of file fit.h.

◆ LTMUnbinnedSig()

template<typename T >
bool o2::math_utils::LTMUnbinnedSig ( const std::vector< T > &  data,
std::vector< size_t > &  index,
std::array< float, 7 > &  params,
float  fracKeepMin,
float  sigTgt,
bool  sorted = false 
)

Compare this function to LTMUnbinned. A target sigma of the distribution can be specified and it will be trimmed to match that target.

Parameters
dataInput vector (unsorted)
indexVector with indices of sorted input data
paramsArray with characteristics of distribution
fracKeepMinMinimum fraction to keep of the input data
sigTgtTarget distribution sigma
sortedFlag if the data is already sorted
Returns
Flag if successfull

Definition at line 648 of file fit.h.

◆ MAD2Sigma()

template<typename T >
T o2::math_utils::MAD2Sigma ( int  np,
T *  y 
)

Definition at line 773 of file fit.h.

◆ median()

template<typename T , typename R = double>
R o2::math_utils::median ( std::vector< T >  v)

median of values in a std::vector

we need to make a copy of the vector since we need to sort it based on this discussion: https://stackoverflow.com/questions/1719070/what-is-the-right-approach-when-using-stl-container-for-median-calculation/1719155#1719155

Todo:
Is there a better way to do this?

Definition at line 519 of file fit.h.

◆ medmadGaus()

template<typename T >
bool o2::math_utils::medmadGaus ( size_t  nBins,
const T *  arr,
const T  xMin,
const T  xMax,
std::array< double, 3 > &  param 
)

fast median estimate of gaussian parameters for histogrammed data

Parameters
[in]nbinssize of the array and number of histogram bins
[in]arrarray with elements
[in]xMinminimum range of the array
[in]xMaxmaximum range of the array
[out]paramreturn paramters of the fit (0-Constant, 1-Mean, 2-Sigma)
Returns
false on failure (empty data)

Definition at line 132 of file fit.h.

◆ Reorder()

template<typename T >
void o2::math_utils::Reorder ( std::vector< T > &  data,
const std::vector< size_t > &  index 
)

Rearranges the input vector in the order given by the index vector

Parameters
dataInput vector
indexIndex vector

Definition at line 625 of file fit.h.

◆ rotateZ() [1/2]

std::tuple< float, float > o2::math_utils::rotateZ ( float  xL,
float  yL,
float  snAlp,
float  csAlp 
)
inline

Definition at line 162 of file Utils.h.

◆ rotateZ() [2/2]

void o2::math_utils::rotateZ ( std::array< float, 3 > &  xy,
float  alpha 
)
inline

Definition at line 172 of file Utils.h.

◆ rotateZd() [1/2]

std::tuple< double, double > o2::math_utils::rotateZd ( double  xL,
double  yL,
double  snAlp,
double  csAlp 
)
inline

Definition at line 167 of file Utils.h.

◆ rotateZd() [2/2]

void o2::math_utils::rotateZd ( std::array< double, 3 > &  xy,
double  alpha 
)
inline

Definition at line 177 of file Utils.h.

◆ rotateZInv()

std::tuple< float, float > o2::math_utils::rotateZInv ( float  xG,
float  yG,
float  snAlp,
float  csAlp 
)
inline

Definition at line 142 of file Utils.h.

◆ rotateZInvd()

std::tuple< double, double > o2::math_utils::rotateZInvd ( double  xG,
double  yG,
double  snAlp,
double  csAlp 
)
inline

Definition at line 147 of file Utils.h.

◆ sector2Angle()

float o2::math_utils::sector2Angle ( int  sect)
inline

Definition at line 193 of file Utils.h.

◆ sector2Angled()

double o2::math_utils::sector2Angled ( int  sect)
inline

Definition at line 198 of file Utils.h.

◆ selKthMin()

template<typename T >
T o2::math_utils::selKthMin ( int  k,
int  np,
T *  arr 
)

Definition at line 718 of file fit.h.

◆ Similarity()

template<class T , unsigned int D1, unsigned int D2, class R >
SMatrix< T, D1, D1, MatRepSym< T, D1 > > o2::math_utils::Similarity ( const SMatrix< T, D1, D2, R > &  lhs,
const SMatrix< T, D2, D2, MatRepSym< T, D2 > > &  rhs 
)
inline

Definition at line 263 of file Cartesian.h.

◆ SortData()

template<typename T >
void o2::math_utils::SortData ( std::vector< T > const &  values,
std::vector< size_t > &  index 
)

Fills the index vector with sorted indices of the input vector. The input vector is not modified (similar to TMath::Sort()).

Parameters
valuesVector to be indexed
indexVector to hold the sorted indices (must have the same size as values)

Definition at line 539 of file fit.h.

◆ to02PiGen()

float o2::math_utils::to02PiGen ( float  phi)
inline

Definition at line 70 of file Utils.h.

◆ to02PiGend()

double o2::math_utils::to02PiGend ( double  phi)
inline

Definition at line 75 of file Utils.h.

◆ toPMPi()

float o2::math_utils::toPMPi ( float  phi)
inline

Definition at line 90 of file Utils.h.

◆ toPMPid()

double o2::math_utils::toPMPid ( double  phi)
inline

Definition at line 95 of file Utils.h.

◆ toPMPiGen()

float o2::math_utils::toPMPiGen ( float  phi)
inline

Definition at line 50 of file Utils.h.

◆ toPMPiGend()

double o2::math_utils::toPMPiGend ( double  phi)
inline

Definition at line 55 of file Utils.h.

Variable Documentation

◆ c

double double & o2::math_utils::c
Initial value:
{
detail::sincos<float>(ang, s, c)
uint32_t c
Definition RawData.h:2

Definition at line 110 of file Utils.h.

◆ csAlp

double double double double double o2::math_utils::csAlp
Initial value:
{
return detail::rotateZ<float>(xL, yL, xG, yG, snAlp, csAlp)

Definition at line 121 of file Utils.h.

◆ s

double & o2::math_utils::s

Definition at line 110 of file Utils.h.

◆ snAlp

double double double double o2::math_utils::snAlp

Definition at line 121 of file Utils.h.

◆ x

double o2::math_utils::x
Initial value:
{
return detail::fastATan2<float>(y, x)
GLint GLenum GLint x
Definition glcorearb.h:403

Definition at line 213 of file Utils.h.

◆ xG

double double & o2::math_utils::xG

Definition at line 121 of file Utils.h.

◆ xL

double double & o2::math_utils::xL

Definition at line 131 of file Utils.h.

◆ y

const double o2::math_utils::y
Initial value:
{
return detail::min<T>(x, y)

Definition at line 224 of file Utils.h.

◆ yG

double o2::math_utils::yG

Definition at line 121 of file Utils.h.

◆ yL

double double double & o2::math_utils::yL

Definition at line 121 of file Utils.h.