Jpp  18.3.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Functions | Variables
JMATH Namespace Reference

Auxiliary classes and methods for mathematical operations. More...

Classes

struct  JAverage< JQuaternion3D >
 Template spacialisation for averaging quaternions. More...
 
struct  JLegendre< JQuaternion3D,(size_t)-1 >
 Template specialisation for function evaluation of Legendre polynome of quaternions for undefined number of degrees. More...
 
struct  JLegendre< JQuaternion3D, N >
 Template specialisation for function evaluation of Legendre polynome of quaternions for defined number of degrees. More...
 
struct  JCalculator
 Auxiliary class for arithmetic operations on objects. More...
 
struct  JGauss_t
 Gauss model. More...
 
struct  JGauss
 Gauss function object. More...
 
struct  JLegendre_t
 Base class for Legendre polynome. More...
 
struct  JLegendre
 Template definition for function evaluation of Legendre polynome. More...
 
struct  JLegendre< JOrdinate_t,(size_t)-1 >
 Template specialisation for function evaluation of of Legendre polynome for undefined number of degrees. More...
 
struct  JLimits
 Auxiliary class for minimum and maximum values for any class. More...
 
struct  JLimits< T, true >
 Template spacialisation of JMATH::JLimits for numerical values. More...
 
struct  JLimits< T, false >
 Template spacialisation of JMATH::JRandom for other data types. More...
 
struct  JMath_t
 Auxiliary class to hide data type specific methods. More...
 
struct  JMath
 Auxiliary base class for aritmetic operations of derived class types. More...
 
struct  JMath< T, JNullType >
 Template base class for data structures with arithmetic capabilities. More...
 
struct  JAverage
 Auxiliary class to determine average of set of values. More...
 
class  JMatrix1D
 1 x 1 matrix More...
 
class  JMatrix1S
 1 x 1 symmetric matrix More...
 
class  JMatrix2D
 2 x 2 matrix More...
 
class  JMatrix2S
 2 x 2 symmetric matrix More...
 
class  JMatrix3D
 3 x 3 matrix More...
 
class  JMatrix3S
 3 x 3 symmetric matrix More...
 
class  JMatrix4D
 4 x 4 matrix More...
 
class  JMatrix4S
 4 x 4 symmetric matrix More...
 
class  JMatrix5D
 5 x 5 matrix More...
 
class  JMatrix5S
 5 x 5 symmetric matrix More...
 
struct  JMatrixND_t
 Basic NxN matrix. More...
 
struct  JMatrixND
 NxN matrix. More...
 
struct  JMatrixNS
 N x N symmetric matrix. More...
 
struct  JModel_t
 Fit model. More...
 
struct  JNumber
 Simple wrapper around template data type to ensure that zero is the default value. More...
 
struct  JPolynome_t
 Polynome model. More...
 
struct  JPolynome
 Polynome function object. More...
 
class  JPower
 Power law function object. More...
 
struct  JRandom
 Template definition of random value generator. More...
 
struct  JRandom< T, true >
 Template spacialisation of JMATH::JRandom for numerical values. More...
 
struct  JRandom< T, false >
 Template spacialisation of JMATH::JRandom for non-numerical data types. More...
 
class  JSVD3D
 Singular value decomposition. More...
 
class  JTrigonometric
 Trigonometric function object for sin and cos. More...
 
struct  JVectorND
 Nx1 matrix. More...
 
struct  JZero
 Auxiliary class to assign zero value. More...
 

Functions

template<class T >
const JCalculator< T, 1 > & operator* (const T &first, const T &second)
 Product evaluation of objects. More...
 
template<class T , int N>
const JCalculator< T, N+1 > & operator* (const T &first, const JCalculator< T, N > &second)
 Recursive product evaluation of objects. More...
 
template<class T , int N>
const JCalculator< T, N+1 > & operator* (const JCalculator< T, N > &first, const T &second)
 Recursive product evaluation of objects. More...
 
template<class T >
T pow (const T &x, const double y)
 Power $ x^{y} $. More...
 
template<class T >
T interpolate (const T &first, const T &second, const double alpha)
 Interpolation between objects. More...
 
template<class T >
std::iterator_traits< T >
::value_type 
getAverage (T __begin, T __end)
 Get average. More...
 
template<class JValue_t , size_t N>
JValue_t getAverage (const JValue_t(&array)[N])
 Get average. More...
 
template<class JElement_t , class JAllocator_t >
JElement_t getAverage (const array_type< JElement_t, JAllocator_t > &buffer)
 Get average. More...
 
template<class T >
std::iterator_traits< T >
::value_type 
getAverage (T __begin, T __end, typename std::iterator_traits< T >::value_type value)
 Get average. More...
 
template<class JValue_t , size_t N>
JValue_t getAverage (const JValue_t(&array)[N], typename JLANG::JClass< JValue_t >::argument_type value)
 Get average. More...
 
template<class JElement_t , class JAllocator_t >
JElement_t getAverage (const array_type< JElement_t, JAllocator_t > &buffer, typename JLANG::JClass< JElement_t >::argument_type value)
 Get average. More...
 
double gauss (const double x, const double sigma)
 Gauss function (normalised to 1 at x = 0). More...
 
double gauss (const double x, const double x0, const double sigma)
 Gauss function (normalised to 1 at x = x0). More...
 
double Gauss (const double x, const double sigma)
 Normalised Gauss function. More...
 
double Gauss (const double x, const double x0, const double sigma)
 Normalised Gauss function. More...
 
double Gamma (const double a, const double x)
 Incomplete gamma function. More...
 
double legendre (const size_t n, const double x)
 Legendre polynome. More...
 
double binomial (const size_t n, const size_t k)
 Binomial function. More...
 
double poisson (const size_t n, const double mu)
 Poisson probability density distribition. More...
 
double Poisson (const size_t n, const double mu)
 Poisson cumulative density distribition. More...
 
void randomize (JMatrix1D *p)
 Randomize matrix. More...
 
void randomize (JMatrix2D *p)
 Randomize matrix. More...
 
void randomize (JMatrix3D *p)
 Randomize matrix. More...
 
void randomize (JMatrix4D *p)
 Randomize matrix. More...
 
void randomize (JMatrix5D *p)
 Randomize matrix. More...
 
void randomize (JMatrix1S *p)
 Randomize matrix. More...
 
void randomize (JMatrix2S *p)
 Randomize matrix. More...
 
void randomize (JMatrix3S *p)
 Randomize matrix. More...
 
void randomize (JMatrix4S *p)
 Randomize matrix. More...
 
void randomize (JMatrix5S *p)
 Randomize matrix. More...
 
long long int factorial (const long long int n)
 Determine factorial. More...
 
long long int factorial (const long long int n, const long long int m)
 Determine combinatorics. More...
 
template<class JFirst_t , class JSecond_t >
bool equals (const JFirst_t &first, const JSecond_t &second, const double precision=std::numeric_limits< double >::min())
 Check equality. More...
 
template<class JFirst_t , class JSecond_t >
double getDistanceSquared (const JFirst_t &first, const JSecond_t &second)
 Get square of distance between objects. More...
 
template<class JFirst_t , class JSecond_t >
double getDistance (const JFirst_t &first, const JSecond_t &second)
 Get distance between objects. More...
 
template<class JFirst_t , class JSecond_t >
double getDot (const JFirst_t &first, const JSecond_t &second)
 Get dot product of objects. More...
 
template<class JFirst_t , class JSecond_t >
double getAngle (const JFirst_t &first, const JSecond_t &second)
 Get space angle between objects. More...
 
template<class JFirst_t , class JSecond_t >
double getPerpDot (const JFirst_t &first, const JSecond_t &second)
 Get perpendicular dot product of objects. More...
 
template<class T >
T getCross (const T &first, const T &second)
 Get cross product of objects. More...
 
template<class T >
std::vector< Tconvolve (const std::vector< T > &input, const std::vector< T > &kernel)
 Convolute data with given kernel. More...
 
template<class T >
T getRandom ()
 Get random value. More...
 
template<class T >
T getRandom (const T min, const T max)
 Get uniformly distributed random value between given limits. More...
 
template<class T >
T getRandom (const T min, const T max, const T precision)
 Get uniformly distributed random value between given limits. More...
 
template<>
double getRandom (const double min, const double max, const double precision)
 Template specialisation for data type double. More...
 
template<class T >
T getZero ()
 Get zero value for a given data type. More...
 
template<>
bool getZero< bool > ()
 Get zero value for bool. More...
 
template<>
float getZero< float > ()
 Get zero value for float. More...
 
template<>
double getZero< double > ()
 Get zero value for double. More...
 
template<>
long double getZero< long double > ()
 Get zero value for long double. More...
 

Variables

static const double PI = acos(-1.0)
 Mathematical constants. More...
 
static const double EULER = 0.577215664901533
 Euler number. More...
 
static const long long int KILOBYTE = 1024
 Computing quantities. More...
 
static const long long int MEGABYTE = KILOBYTE*KILOBYTE
 Number of bytes in a kilo-byte. More...
 
static const long long int GIGABYTE = MEGABYTE*KILOBYTE
 Number of bytes in a mega-byte. More...
 
static const JZero zero
 Function object to assign zero value. More...
 

Detailed Description

Auxiliary classes and methods for mathematical operations.

Author
mdejong

Function Documentation

template<class T >
const JCalculator<T, 1>& JMATH::operator* ( const T first,
const T second 
)
inline

Product evaluation of objects.

Multiply objects.

Parameters
firstfirst object
secondsecond object
Returns
calculator

Definition at line 53 of file JCalculator.hh.

54  {
56 
58  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
Auxiliary class for arithmetic operations on objects.
Definition: JCalculator.hh:18
template<class T , int N>
const JCalculator<T, N+1>& JMATH::operator* ( const T first,
const JCalculator< T, N > &  second 
)
inline

Recursive product evaluation of objects.

Parameters
firstfirst object
secondsecond object
Returns
calculator

Definition at line 69 of file JCalculator.hh.

70  {
72 
74  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
Auxiliary class for arithmetic operations on objects.
Definition: JCalculator.hh:18
template<class T , int N>
const JCalculator<T, N+1>& JMATH::operator* ( const JCalculator< T, N > &  first,
const T second 
)
inline

Recursive product evaluation of objects.

Parameters
firstfirst object
secondsecond object
Returns
calculator

Definition at line 85 of file JCalculator.hh.

86  {
87  JCalculator<T, N+1>::calculator.mul(first, second);
88 
90  }
Auxiliary class for arithmetic operations on objects.
Definition: JCalculator.hh:18
template<class T >
T JMATH::pow ( const T x,
const double  y 
)
inline

Power $ x^{y} $.

Parameters
xvalue
ypower
Returns
result

Definition at line 97 of file JMath.hh.

98  {
99  using namespace JPP;
100 
101  return JMath_t::pow(x, y, JBool<JClass<T>::is_primitive>());
102  }
static T pow(const T &x, const double y, const JLANG::JBool< true > option)
Power .
Definition: JMath.hh:65
template<class T >
T JMATH::interpolate ( const T first,
const T second,
const double  alpha 
)
inline

Interpolation between objects.

The result is equal to result = (1 - alpha) * (first) + (alpha) * (second).

Parameters
firstfirst object
secondsecond object
alphainterpolation factor [0, 1]
Returns
result

Definition at line 397 of file JMath.hh.

400  {
401  return T(first).interpolate(second, alpha);
402  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
do set_variable OUTPUT_DIRECTORY $WORKDIR T
template<class T >
std::iterator_traits<T>::value_type JMATH::getAverage ( T  __begin,
T  __end 
)

Get average.

Parameters
__beginbegin of data
__endend of data
Returns
average value

Definition at line 494 of file JMath.hh.

495  {
496  typedef typename std::iterator_traits<T>::value_type value_type;
497 
498  return JAverage<value_type>(__begin, __end);
499  }
Auxiliary class to determine average of set of values.
Definition: JMath.hh:409
template<class JValue_t , size_t N>
JValue_t JMATH::getAverage ( const JValue_t(&)  array[N])
inline

Get average.

Parameters
arrayc-array of values
Returns
average value

Definition at line 509 of file JMath.hh.

510  {
511  typedef JValue_t value_type;
512 
513  return JAverage<value_type>((const value_type*) array, (const value_type*) array + N);
514  }
then usage $script< input file >[option[primary[working directory]]] nWhere option can be N
Definition: JMuonPostfit.sh:40
Auxiliary class to determine average of set of values.
Definition: JMath.hh:409
template<class JElement_t , class JAllocator_t >
JElement_t JMATH::getAverage ( const array_type< JElement_t, JAllocator_t > &  buffer)

Get average.

Parameters
bufferinput data
Returns
average value

Definition at line 524 of file JMath.hh.

525  {
526  return JAverage<JElement_t>(buffer.begin(), buffer.end());
527  }
Auxiliary class to determine average of set of values.
Definition: JMath.hh:409
template<class T >
std::iterator_traits<T>::value_type JMATH::getAverage ( T  __begin,
T  __end,
typename std::iterator_traits< T >::value_type  value 
)

Get average.

Parameters
__beginbegin of data
__endend of data
valuedefault value
Returns
average value

Definition at line 539 of file JMath.hh.

540  {
541  try {
542  return getAverage(__begin, __end);
543  }
544  catch(const std::exception&) {
545  return value;
546  }
547  }
std::iterator_traits< T >::value_type getAverage(T __begin, T __end)
Get average.
Definition: JMath.hh:494
template<class JValue_t , size_t N>
JValue_t JMATH::getAverage ( const JValue_t(&)  array[N],
typename JLANG::JClass< JValue_t >::argument_type  value 
)
inline

Get average.

Parameters
arrayc-array of values
valuedefault value
Returns
average value

Definition at line 558 of file JMath.hh.

559  {
560  try {
561  return getAverage(array);
562  }
563  catch(const std::exception&) {
564  return value;
565  }
566  }
std::iterator_traits< T >::value_type getAverage(T __begin, T __end)
Get average.
Definition: JMath.hh:494
template<class JElement_t , class JAllocator_t >
JElement_t JMATH::getAverage ( const array_type< JElement_t, JAllocator_t > &  buffer,
typename JLANG::JClass< JElement_t >::argument_type  value 
)

Get average.

Parameters
bufferinput data
valuedefault value
Returns
average value

Definition at line 577 of file JMath.hh.

578  {
579  try {
580  return getAverage(buffer);
581  }
582  catch(const std::exception&) {
583  return value;
584  }
585  }
std::iterator_traits< T >::value_type getAverage(T __begin, T __end)
Get average.
Definition: JMath.hh:494
double JMATH::gauss ( const double  x,
const double  sigma 
)
inline

Gauss function (normalised to 1 at x = 0).

Parameters
xx
sigmasigma
Returns
function value

Definition at line 33 of file JMath/JMathSupportkit.hh.

34  {
35  const double u = x / sigma;
36 
37  if (fabs(u) < 20.0)
38  return exp(-0.5*u*u);
39  else
40  return 0.0;
41  }
const double sigma[]
Definition: JQuadrature.cc:74
double u[N+1]
Definition: JPolint.hh:865
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
double JMATH::gauss ( const double  x,
const double  x0,
const double  sigma 
)
inline

Gauss function (normalised to 1 at x = x0).

Parameters
xx
x0central value
sigmasigma
Returns
function value

Definition at line 52 of file JMath/JMathSupportkit.hh.

53  {
54  return gauss(x - x0, sigma);
55  }
const double sigma[]
Definition: JQuadrature.cc:74
then set_variable NUMBER_OF_TESTS else set_variable NUMBER_OF_TESTS fi function gauss()
double JMATH::Gauss ( const double  x,
const double  sigma 
)
inline

Normalised Gauss function.

Parameters
xx
sigmasigma
Returns
function value

Definition at line 65 of file JMath/JMathSupportkit.hh.

66  {
67  return gauss(x, sigma) / sqrt(2.0*PI) / sigma;
68  }
const double sigma[]
Definition: JQuadrature.cc:74
static const double PI
Mathematical constants.
then set_variable NUMBER_OF_TESTS else set_variable NUMBER_OF_TESTS fi function gauss()
double JMATH::Gauss ( const double  x,
const double  x0,
const double  sigma 
)
inline

Normalised Gauss function.

Parameters
xx
x0central value
sigmasigma
Returns
function value

Definition at line 79 of file JMath/JMathSupportkit.hh.

80  {
81  return Gauss(x - x0, sigma);
82  }
const double sigma[]
Definition: JQuadrature.cc:74
double Gauss(const double x, const double sigma)
Normalised Gauss function.
double JMATH::Gamma ( const double  a,
const double  x 
)
inline

Incomplete gamma function.

Source code is taken from reference: Numerical Recipes in C++, W.H. Press, S.A. Teukolsky, W.T. Vetterling and B.P. Flannery, Cambridge University Press.

Parameters
aa
xx
Returns
function value

Definition at line 96 of file JMath/JMathSupportkit.hh.

97  {
98  using namespace std;
99 
100  const int max = 1000000;
101 
102  if (x < 0.0) { THROW(JValueOutOfRange, "x < 0 " << x); }
103  if (a <= 0.0) { THROW(JValueOutOfRange, "a <= 0 " << a); }
104 
105  const double gln = lgamma(a);
106 
107  if (x < a + 1.0) {
108 
109  if (x < 0.0) {
110  THROW(JValueOutOfRange, "x <= 0 " << x);
111  }
112 
113  if (x == 0.0) {
114  return 0.0;
115  }
116 
117  double ap = a;
118  double sum = 1.0 / a;
119  double del = sum;
120 
121  for (int i = 0; i != max; ++i) {
122 
123  ap += 1.0;
124  del *= x/ap;
125  sum += del;
126 
127  if (fabs(del) < fabs(sum)*numeric_limits<double>::epsilon()) {
128  return sum*exp(-x + a*log(x) - gln);
129  }
130  }
131 
132  THROW(JValueOutOfRange, "i == " << max);
133 
134  } else {
135 
136  const double FPMIN = numeric_limits<double>::min() / numeric_limits<double>::epsilon();
137 
138  double b = x + 1.0 - a;
139  double c = 1.0 / FPMIN;
140  double d = 1.0 / b;
141  double h = d;
142 
143  for (int i = 1; i != max; ++i) {
144 
145  const double an = -i * (i-a);
146 
147  b += 2.0;
148  d = an*d + b;
149 
150  if (fabs(d) < FPMIN) {
151  d = FPMIN;
152  }
153 
154  c = b + an/c;
155 
156  if (fabs(c) < FPMIN) {
157  c = FPMIN;
158  }
159 
160  d = 1.0/d;
161 
162  const double del = d*c;
163 
164  h *= del;
165 
166  if (fabs(del - 1.0) < numeric_limits<double>::epsilon()) {
167  return 1.0 - exp(-x + a*log(x) - gln) * h;
168  }
169  }
170 
171  THROW(JValueOutOfRange, "i == " << max);
172  }
173  }
then set_variable DIR else fatal Wrong number of arguments fi for INPUT_FILE in ls rt $DIR stage * log
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:712
then JCalibrateToT a
Definition: JTuneHV.sh:113
$WORKDIR ev_configure_dqsimulator txt echo process $DQ_SIMULATOR $i $SOURCE_HOST[$index] csh c(setenv ROOTSYS $ROOTSYS &&source $JPP_DIR/setenv.csh $JPP_DIR &&($DQ_SIMULATOR\-u\$NAME\$\-H\$SERVER\$\-M\$LOGGER\$\-d $DEBUG</dev/null > &/dev/null &))'
then JMuonMCEvt f $INPUT_FILE o $INTERMEDIATE_FILE d
Definition: JMuonPath.sh:47
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
const double epsilon
Definition: JQuadrature.cc:21
double JMATH::legendre ( const size_t  n,
const double  x 
)
inline

Legendre polynome.

Parameters
ndegree
xx
Returns
function value

Definition at line 183 of file JMath/JMathSupportkit.hh.

184  {
185  switch (n) {
186 
187  case 0:
188  return 1.0;
189 
190  case 1:
191  return x;
192 
193  default:
194  {
195  double p0;
196  double p1 = 1.0;
197  double p2 = x;
198 
199  for (size_t i = 2; i <= n; ++i) {
200  p0 = p1;
201  p1 = p2;
202  p2 = ((2*i-1) * x*p1 - (i-1) * p0) / i;
203  }
204 
205  return p2;
206  }
207  }
208  }
TPaveText * p1
const int n
Definition: JPolint.hh:786
p2
Definition: module-Z:fit.sh:74
double JMATH::binomial ( const size_t  n,
const size_t  k 
)
inline

Binomial function.

Parameters
nn
kk
Returns
function value

Definition at line 218 of file JMath/JMathSupportkit.hh.

219  {
220  if (n == 0 || n < k) {
221  return 0.0;
222  }
223 
224  if (k == 0 || n == k) {
225  return 1.0;
226  }
227 
228  const int k1 = std::min(k, n - k);
229  const int k2 = n - k1;
230 
231  double value = k2 + 1;
232 
233  for (int i = k1; i != 1; --i) {
234  value *= (double) (k2 + i) / (double) i;
235  }
236 
237  return value;
238  }
then fatal No hydrophone data file $HYDROPHONE_TXT fi sort gr k
const int n
Definition: JPolint.hh:786
double JMATH::poisson ( const size_t  n,
const double  mu 
)
inline

Poisson probability density distribition.

Parameters
nnumber of occurences
muexpectation value
Returns
probability

Definition at line 248 of file JMath/JMathSupportkit.hh.

249  {
250  using namespace std;
251 
252  if (mu > 0.0) {
253 
254  if (n > 0)
255  return exp(n*log(mu) - lgamma(n+1) - mu);
256  else
257  return exp(-mu);
258  } else if (mu == 0.0) {
259 
260  return (n == 0 ? 1.0 : 0.0);
261  }
262 
263  THROW(JValueOutOfRange, "mu <= 0 " << mu);
264  }
then set_variable DIR else fatal Wrong number of arguments fi for INPUT_FILE in ls rt $DIR stage * log
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:712
const int n
Definition: JPolint.hh:786
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
double JMATH::Poisson ( const size_t  n,
const double  mu 
)
inline

Poisson cumulative density distribition.

Parameters
nnumber of occurences
muexpectation value
Returns
probability

Definition at line 274 of file JMath/JMathSupportkit.hh.

275  {
276  return 1.0 - Gamma(n + 1, mu);
277  }
double Gamma(const double a, const double x)
Incomplete gamma function.
const int n
Definition: JPolint.hh:786
void JMATH::randomize ( JMatrix1D p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 33 of file JMathTestkit.hh.

34  {
35  new (p) JMatrix1D(getRandom<double>(-1.0, +1.0));
36  }
1 x 1 matrix
Definition: JMatrix1D.hh:32
void JMATH::randomize ( JMatrix2D p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 44 of file JMathTestkit.hh.

45  {
46  new (p) JMatrix2D(getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
47  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0));
48  }
2 x 2 matrix
void JMATH::randomize ( JMatrix3D p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 56 of file JMathTestkit.hh.

57  {
58  new (p) JMatrix3D(getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
59  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
60  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0));
61  }
3 x 3 matrix
void JMATH::randomize ( JMatrix4D p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 69 of file JMathTestkit.hh.

70  {
71  new (p) JMatrix4D(getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
72  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
73  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
74  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0));
75  }
4 x 4 matrix
Definition: JMatrix4D.hh:33
void JMATH::randomize ( JMatrix5D p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 83 of file JMathTestkit.hh.

84  {
85  new (p) JMatrix5D(getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
86  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
87  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
88  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
89  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0));
90  }
5 x 5 matrix
Definition: JMatrix5D.hh:33
void JMATH::randomize ( JMatrix1S p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 98 of file JMathTestkit.hh.

99  {
100  new (p) JMatrix1S(getRandom<double>(-1.0, +1.0));
101  }
1 x 1 symmetric matrix
Definition: JMatrix1S.hh:26
void JMATH::randomize ( JMatrix2S p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 109 of file JMathTestkit.hh.

110  {
111  new (p) JMatrix2S(getRandom<double>(-1.0,+1.0),
112  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0));
113 
114  p->a01 = p->a10;
115  }
2 x 2 symmetric matrix
Definition: JMatrix2S.hh:26
void JMATH::randomize ( JMatrix3S p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 123 of file JMathTestkit.hh.

124  {
125  new (p) JMatrix3S(getRandom<double>(-1.0,+1.0),
126  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
127  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0));
128 
129  p->a01 = p->a10;
130  p->a02 = p->a20;
131  p->a12 = p->a21;
132  }
3 x 3 symmetric matrix
Definition: JMatrix3S.hh:29
void JMATH::randomize ( JMatrix4S p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 140 of file JMathTestkit.hh.

141  {
142  new (p) JMatrix4S(getRandom<double>(-1.0,+1.0),
143  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
144  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
145  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0));
146 
147  p->a01 = p->a10;
148  p->a02 = p->a20;
149  p->a03 = p->a30;
150  p->a12 = p->a21;
151  p->a13 = p->a31;
152  p->a23 = p->a32;
153  }
4 x 4 symmetric matrix
Definition: JMatrix4S.hh:26
void JMATH::randomize ( JMatrix5S p)
inline

Randomize matrix.

Parameters
ppointer to valid object

Definition at line 161 of file JMathTestkit.hh.

162  {
163  new (p) JMatrix5S(getRandom<double>(-1.0,+1.0),
164  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
165  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
166  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0),
167  getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0), getRandom<double>(-1.0,+1.0));
168 
169  p->a01 = p->a10;
170  p->a02 = p->a20;
171  p->a03 = p->a30;
172  p->a04 = p->a40;
173  p->a12 = p->a21;
174  p->a13 = p->a31;
175  p->a14 = p->a41;
176  p->a23 = p->a32;
177  p->a24 = p->a42;
178  p->a34 = p->a43;
179  }
5 x 5 symmetric matrix
Definition: JMatrix5S.hh:26
long long int JMATH::factorial ( const long long int  n)
inline

Determine factorial.

Parameters
nnumber
Returns
factorial (= n!)

Definition at line 42 of file JMathToolkit.hh.

43  {
44  if (n < 0) {
45  THROW(JValueOutOfRange, "JMATH::factorial(): invalid argument " << n);
46  }
47 
48  if (n != 1)
49  return n * factorial(n - 1);
50  else
51  return 1;
52  }
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:712
long long int factorial(const long long int n)
Determine factorial.
Definition: JMathToolkit.hh:42
const int n
Definition: JPolint.hh:786
long long int JMATH::factorial ( const long long int  n,
const long long int  m 
)
inline

Determine combinatorics.

Parameters
nnumber
mnumber
Returns
n!/(m!*(n-m)!)

Definition at line 62 of file JMathToolkit.hh.

63  {
64  if (n < 0 || m < 0 || n < m) {
65  THROW(JValueOutOfRange, "JMATH::factorial(): invalid argument " << n << ' ' << m);
66  }
67 
68  if (n == m)
69  return 1;
70  else if (m == 0)
71  return 1;
72  else
73  return factorial(n - 1, m - 1) + factorial(n - 1, m);
74  }
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:712
long long int factorial(const long long int n)
Determine factorial.
Definition: JMathToolkit.hh:42
const int n
Definition: JPolint.hh:786
template<class JFirst_t , class JSecond_t >
bool JMATH::equals ( const JFirst_t &  first,
const JSecond_t &  second,
const double  precision = std::numeric_limits<double>::min() 
)
inline

Check equality.

Parameters
firstfirst object
secondsecond object
precisionprecision
Returns
true if two objects are equals; else false

Definition at line 86 of file JMathToolkit.hh.

89  {
90  return first.equals(second, precision);
91  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
template<class JFirst_t , class JSecond_t >
double JMATH::getDistanceSquared ( const JFirst_t &  first,
const JSecond_t &  second 
)
inline

Get square of distance between objects.

Parameters
firstfirst object
secondsecond object
Returns
square of distance

Definition at line 102 of file JMathToolkit.hh.

104  {
105  return first.getDistanceSquared(second);
106  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
template<class JFirst_t , class JSecond_t >
double JMATH::getDistance ( const JFirst_t &  first,
const JSecond_t &  second 
)
inline

Get distance between objects.

Parameters
firstfirst object
secondsecond object
Returns
distance

Definition at line 117 of file JMathToolkit.hh.

119  {
120  return first.getDistance(second);
121  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
template<class JFirst_t , class JSecond_t >
double JMATH::getDot ( const JFirst_t &  first,
const JSecond_t &  second 
)
inline

Get dot product of objects.

Parameters
firstfirst object
secondsecond object
Returns
dot product

Definition at line 132 of file JMathToolkit.hh.

134  {
135  return first.getDot(second);
136  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
template<class JFirst_t , class JSecond_t >
double JMATH::getAngle ( const JFirst_t &  first,
const JSecond_t &  second 
)
inline

Get space angle between objects.

Parameters
firstfirst object
secondsecond object
Returns
angle [deg]

Definition at line 147 of file JMathToolkit.hh.

149  {
150  const double dot = getDot(first,second);
151 
152  if (dot >= +1.0)
153  return 0.0;
154  else if (dot <= -1.0)
155  return 180.0;
156  else
157  return acos(dot) * 180.0 / acos(-1.0);
158  }
double getDot(const JNeutrinoDirection &first, const JNeutrinoDirection &second)
Dot product.
Definition: JAstronomy.hh:674
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
template<class JFirst_t , class JSecond_t >
double JMATH::getPerpDot ( const JFirst_t &  first,
const JSecond_t &  second 
)
inline

Get perpendicular dot product of objects.

Parameters
firstfirst object
secondsecond object
Returns
perpendicular dot product

Definition at line 169 of file JMathToolkit.hh.

171  {
172  return first.getPerpDot(second);
173  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
template<class T >
T JMATH::getCross ( const T first,
const T second 
)
inline

Get cross product of objects.

Parameters
firstfirst object
secondsecond object
Returns
cross product

Definition at line 184 of file JMathToolkit.hh.

186  {
187  return T().getCross(first, second);
188  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
do set_variable OUTPUT_DIRECTORY $WORKDIR T
template<class T >
std::vector<T> JMATH::convolve ( const std::vector< T > &  input,
const std::vector< T > &  kernel 
)
inline

Convolute data with given kernel.

Parameters
inputinput data
kernelconvolution kernel
Returns
convoluted data

Definition at line 199 of file JMathToolkit.hh.

200  {
201  const size_t k = kernel.size();
202  const size_t n = input .size() - k + 1;
203 
204  std::vector<T> out(n, getZero<T>());
205 
206  for (size_t i = 0; i != n; ++i) {
207  for (size_t j = 0; j != k; ++j) {
208  out[i] += input[i + j] * kernel[k - j - 1];
209  }
210  }
211 
212  return out;
213  }
then fatal No hydrophone data file $HYDROPHONE_TXT fi sort gr k
const int n
Definition: JPolint.hh:786
int j
Definition: JPolint.hh:792
template<class T >
T JMATH::getRandom ( )
inline

Get random value.

Returns
random value

Definition at line 113 of file JRandom.hh.

114  {
115  return JRandom<T>::getRandom();
116  }
T getRandom()
Get random value.
Definition: JRandom.hh:113
template<class T >
T JMATH::getRandom ( const T  min,
const T  max 
)
inline

Get uniformly distributed random value between given limits.

Parameters
minminimal value
maxmaximal value
Returns
random value

Definition at line 127 of file JRandom.hh.

129  {
130  return JRandom<T, true>::getRandom(min, max);
131  }
T getRandom()
Get random value.
Definition: JRandom.hh:113
template<class T >
T JMATH::getRandom ( const T  min,
const T  max,
const T  precision 
)
inline

Get uniformly distributed random value between given limits.

Parameters
minminimal value
maxmaximal value
precisionprecision
Returns
random value

Definition at line 143 of file JRandom.hh.

146  {
147  return JRandom<T, true>::getRandom(min, max);
148  }
T getRandom()
Get random value.
Definition: JRandom.hh:113
template<>
double JMATH::getRandom ( const double  min,
const double  max,
const double  precision 
)
inline

Template specialisation for data type double.

Definition at line 155 of file JRandom.hh.

158  {
159  double value = getRandom<double>(min, max);
160 
161  if (precision != 0.0) {
162 
163  const long long int ip = (long long int) (1.0 / precision);
164 
165  value = ((double) ((long long int) (value * ip))) / ip;
166  }
167 
168  return value;
169  }
template<class T >
T JMATH::getZero ( )
inline

Get zero value for a given data type.

The default implementation of this method returns an object which is created with the default constructor. This method should be specialised if this value does not correspond to the equivalent of a zero result.

Returns
zero

Definition at line 26 of file JZero.hh.

27  {
28  return T();
29  }
do set_variable OUTPUT_DIRECTORY $WORKDIR T
template<>
bool JMATH::getZero< bool > ( )
inline

Get zero value for bool.

Returns
false

Definition at line 37 of file JZero.hh.

38  {
39  return false;
40  }
template<>
float JMATH::getZero< float > ( )
inline

Get zero value for float.

Returns
zero

Definition at line 48 of file JZero.hh.

49  {
50  return float(0.0);
51  }
template<>
double JMATH::getZero< double > ( )
inline

Get zero value for double.

Returns
zero

Definition at line 60 of file JZero.hh.

61  {
62  return double(0.0);
63  }
template<>
long double JMATH::getZero< long double > ( )
inline

Get zero value for long double.

Returns
zero

Definition at line 72 of file JZero.hh.

73  {
74  return (long double)(0.0);
75  }

Variable Documentation

const double JMATH::PI = acos(-1.0)
static

Mathematical constants.

pi

Definition at line 20 of file JMath/JConstants.hh.

const double JMATH::EULER = 0.577215664901533
static

Euler number.

Definition at line 21 of file JMath/JConstants.hh.

const long long int JMATH::KILOBYTE = 1024
static

Computing quantities.

Definition at line 26 of file JMath/JConstants.hh.

const long long int JMATH::MEGABYTE = KILOBYTE*KILOBYTE
static

Number of bytes in a kilo-byte.

Definition at line 27 of file JMath/JConstants.hh.

const long long int JMATH::GIGABYTE = MEGABYTE*KILOBYTE
static

Number of bytes in a mega-byte.

Definition at line 28 of file JMath/JConstants.hh.

const JZero JMATH::zero
static

Function object to assign zero value.

Definition at line 105 of file JZero.hh.