Theoretica
A C++ numerical and automatic mathematical library
real_analysis.h
Go to the documentation of this file.
1 
5 
6 #ifndef THEORETICA_COMMON_H
7 #define THEORETICA_COMMON_H
8 
9 #include "./core_traits.h"
10 #include "./constants.h"
11 #include "./error.h"
12 
13 
14 namespace theoretica {
15 
16 
19  return x;
20  }
21 
22 
24  template<typename Type, typename = std::enable_if<is_real_type<Type>::value>>
25  TH_CONSTEXPR inline Type conjugate(Type x) {
26  return x;
27  }
28 
29 
36  return x * x;
37  }
38 
39 
45  TH_CONSTEXPR inline real cube(real x) {
46  return x * x * x;
47  }
48 
49 
54  template<typename UnsignedIntType = uint64_t>
55  inline UnsignedIntType isqrt(UnsignedIntType n) {
56 
57  // Upper bound
58  UnsignedIntType upper = n + 1;
59 
60  // Lower bound
61  UnsignedIntType lower = 0;
62 
63  // Carry for safe long division
64  UnsignedIntType c = 0;
65 
66  while(lower != upper - 1) {
67 
68  // Compute carry for long division by 2
69  c = ((lower % 2 != 0) && (upper % 2 != 0)) ? 1 : 0;
70 
71  // Safer division by 2 for big numbers
72  const UnsignedIntType m = (lower >> 1) + (upper >> 1) + c;
73 
74  // Using division instead of multiplication avoids
75  // overflows which would remove significant bits
76  const UnsignedIntType q = n / m;
77 
78  if(m > q)
79  upper = m;
80  else if(m < q)
81  lower = m;
82  else
83  return m;
84  }
85 
86  return lower;
87  }
88 
89 
94  template<typename UnsignedIntType = uint64_t>
95  inline UnsignedIntType icbrt(UnsignedIntType n) {
96 
97  // Upper bound
98  UnsignedIntType upper = n + 1;
99 
100  // Lower bound
101  UnsignedIntType lower = 0;
102 
103  // Carry for safe long division
104  UnsignedIntType c = 0;
105 
106  while(lower != upper - 1) {
107 
108  // Compute carry for long division by 2
109  c = ((lower % 2 != 0) && (upper % 2 != 0)) ? 1 : 0;
110 
111  // Safer division by 2 for big numbers
112  const UnsignedIntType m = (lower >> 1) + (upper >> 1) + c;
113 
114  // Using division instead of multiplication avoids
115  // overflows which would remove significant bits
116  const UnsignedIntType q = (n / m) / m;
117 
118  if(m > q)
119  upper = m;
120  else if(m < q)
121  lower = m;
122  else
123  return m;
124  }
125 
126  return lower;
127  }
128 
129 
142  inline real sqrt(real x) {
143 
144  if(x < 0) {
145  TH_MATH_ERROR("sqrt", x, OUT_OF_DOMAIN);
146  return nan();
147  }
148 
149 #ifdef THEORETICA_X86
150 
151  #ifdef MSVC_ASM
152  __asm {
153  fld x
154  fsqrt
155  }
156  #else
157  asm("fsqrt" : "+t"(x));
158  #endif
159 
160  return x;
161 #else
162 
163  if(x < 1) {
164 
165  if(x == 0)
166  return 0;
167 
168  // Approximate sqrt(x) between 0 and 1
169  // The root of the inverse is the inverse of the root
170  // !!! Possible precision problems with smaller numbers
171  return 1.0 / sqrt(1.0 / x);
172  }
173 
174  // Approximate sqrt(x) using Newton-Raphson
175  real y = x;
176  unsigned int i = 0;
177 
178  while((square(y) - x) > OPTIMIZATION_TOL && i < OPTIMIZATION_NEWTON_ITER) {
179  y = (y + x / y) / 2.0;
180  i++;
181  }
182 
183  return y;
184 #endif
185 
186  }
187 
188 
199  inline real cbrt(real x) {
200 
201  if(x < 1) {
202 
203  if(x == 0)
204  return 0;
205 
206  // cbrt(x) is odd
207  if(x < 0)
208  return -cbrt(-x);
209 
210  // Approximate cbrt between 0 and 1
211  // The root of the inverse is the inverse of the root
212  // !!! Possible precision problems with smaller numbers
213  return 1.0 / cbrt(1.0 / x);
214  }
215 
216  // Approximate cbrt(x) using Newton-Raphson
217  real y = x;
218  unsigned int i = 0;
219 
220  while((cube(y) - x) > OPTIMIZATION_TOL && i < OPTIMIZATION_NEWTON_ITER) {
221  y = (y * 2.0 + x / (y * y)) / 3.0;
222  i++;
223  }
224 
225  return y;
226  }
227 
228 
235  inline real abs(real x) {
236 
237 #ifdef THEORETICA_X86
238 
239  #ifdef MSVC_ASM
240  __asm {
241  fld x
242  fabs
243  }
244  #else
245  asm("fabs" : "+t"(x));
246  #endif
247 
248  return x;
249 #else
250  return x >= 0 ? x : -x;
251 #endif
252 
253  }
254 
255 
259  inline int sgn(real x) {
260  return (x > 0) ? 1 : (x < 0 ? -1 : 0);
261  }
262 
263 
271  TH_CONSTEXPR inline int floor(real x) {
272 
273  if(x < 0 && x > -1)
274  return -1;
275 
276  // Compute the biggest integer number
277  // that is smaller than x
278  return x - (int(x) % 1);
279  }
280 
281 
288  inline real fract(real x) {
289  return abs(x - floor(x));
290  }
291 
292 
300  inline real max(real x, real y) {
301 
302  #ifdef THEORETICA_BRANCHLESS
303  return (x + y + abs(x - y)) / 2.0;
304  #else
305  return x > y ? x : y;
306  #endif
307  }
308 
315  template<typename T>
316  inline T max(T x, T y) {
317  return x > y ? x : y;
318  }
319 
320 
328  inline real min(real x, real y) {
329 
330  #ifdef THEORETICA_BRANCHLESS
331  return (x + y - abs(x - y)) / 2.0;
332  #else
333  return x > y ? y : x;
334  #endif
335  }
336 
343  template<typename T>
344  inline T min(T x, T y) {
345  return x > y ? y : x;
346  }
347 
348 
355  inline real clamp(real x, real a, real b) {
356  return x > b ? b : (x < a ? a : x);
357  }
358 
359 
368  template<typename T>
369  inline T clamp(T x, T a, T b) {
370  return x > b ? b : (x < a ? a : x);
371  }
372 
373 
374  // x86 instruction wrappers
375 
376 #ifdef THEORETICA_X86
377 
379  inline real fyl2x(real x, real y) {
380 
381  #ifdef MSVC_ASM
382 
383  __asm {
384  fld y
385  fld x
386  fyl2x
387  }
388 
389  #else
390  double r;
391  asm ("fyl2x" : "=t"(r) : "0"(x), "u"(y) : "st(1)");
392  return r;
393  #endif
394  }
395 
396 
401  inline real f2xm1(real x) {
402 
403  #ifdef MSVC_ASM
404 
405  __asm {
406  fld x
407  f2xm1
408  }
409 
410  #else
411  asm("f2xm1" : "+t"(x));
412  return x;
413  #endif
414  }
415 
416 #endif
417 
425  inline real log2(real x) {
426 
427  if(x <= 0) {
428 
429  if(x == 0) {
430  TH_MATH_ERROR("log2", x, OUT_OF_RANGE);
431  return -inf();
432  }
433 
434  TH_MATH_ERROR("log2", x, OUT_OF_DOMAIN);
435  return nan();
436  }
437 
438 #ifdef THEORETICA_X86
439 
440  // Approximate the binary logarithm of x by
441  // exploiting x86 Assembly instructions
442  return fyl2x(x, 1.0);
443 #else
444 
445  // Domain reduction to [1, +inf)
446  if(x < 1)
447  return -log2(1 / x);
448 
449  // Compute the biggest power of 2
450  // so that x <= 2^i
451  int i = 0;
452  while(x > (1 << i))
453  i++;
454 
455  // Domain reduction to [1, 2]
456  x /= (1 << i);
457 
458  // Use the Taylor expansion of the logarithm
459  // ln(1 + x) = \sum_k^n (-1)^(k+1) x^k / k
460  real log_z = 0;
461  const real z = x - 1;
462 
463  // Exact powers of 2 don't need further computation
464  if(abs(z) > MACH_EPSILON) {
465 
466  int sign = 1;
467  real pow_z = z;
468  log_z = z;
469 
470  for (int i = 2; i <= 24; ++i) {
471 
472  pow_z *= z;
473  sign *= -1;
474 
475  log_z += sign * pow_z / i;
476  }
477  }
478 
479  return i + (log_z / LN2);
480 #endif
481  }
482 
483 
491  inline real log10(real x) {
492 
493  if(x <= 0) {
494 
495  if(x == 0) {
496  TH_MATH_ERROR("log10", x, OUT_OF_RANGE);
497  return -inf();
498  }
499 
500  TH_MATH_ERROR("log10", x, OUT_OF_DOMAIN);
501  return nan();
502  }
503 
504 #ifdef THEORETICA_X86
505 
506  // Approximate the binary logarithm of x by
507  // exploiting x86 Assembly instructions
508  return fyl2x(x, 1.0 / LOG210);
509 #else
510  return log2(x) / LOG210;
511 #endif
512  }
513 
514 
522  inline real ln(real x) {
523 
524  if(x <= 0) {
525 
526  if(x == 0) {
527  TH_MATH_ERROR("ln", x, OUT_OF_RANGE);
528  return -inf();
529  }
530 
531  TH_MATH_ERROR("ln", x, OUT_OF_DOMAIN);
532  return nan();
533  }
534 
535 #ifdef THEORETICA_X86
536 
537  // Approximate the binary logarithm of x by
538  // exploiting x86 Assembly instructions
539  return fyl2x(x, 1.0 / LOG2E);
540 #else
541  return log2(x) / LOG2E;
542 #endif
543  }
544 
545 
550  template<typename UnsignedIntType = uint64_t>
551  inline UnsignedIntType ilog2(UnsignedIntType x) {
552 
553  if(x == 0) {
554  TH_MATH_ERROR("ilog2", x, OUT_OF_RANGE);
555  return -inf();
556  }
557 
558  UnsignedIntType bit = 0;
559 
560  // Find the highest set bit
561  for (int i = (sizeof(UnsignedIntType) * 8 - 1); i > 0; --i) {
562  if(x & ((UnsignedIntType) 1 << i)) {
563  bit = i;
564  break;
565  }
566  }
567 
568  return bit;
569  }
570 
571 
577  template<typename UnsignedIntType = uint64_t>
578  inline UnsignedIntType pad2(UnsignedIntType x) {
579 
580  UnsignedIntType bit = 0;
581 
582  for (int i = (sizeof(UnsignedIntType) * 8 - 1); i > 0; --i) {
583  if(x & ((UnsignedIntType) 1 << i)) {
584 
585  // Find the highest set bit
586  bit = i;
587 
588  // Check if x is a power of 2
589  for (int j = 0; j < i; ++j)
590  if(x & ((UnsignedIntType) 1 << j))
591  return (1 << (bit + 1));
592  }
593  }
594 
595  return (1 << bit);
596  }
597 
598 
603  template<typename T = real>
604  TH_CONSTEXPR inline T pow(T x, int n) {
605 
606  if(n > 0) {
607 
608  T res = x;
609  T x_sqr = x * x;
610  int i = 1;
611 
612  // Self-multiply up to biggest power of 2
613  for (; i < (n / 2); i *= 2)
614  res = res * res;
615 
616  // Multiply by x^2 for remaining even powers
617  for (; i < (n - 1); i += 2)
618  res = res * x_sqr;
619 
620  // Multiply for remaining powers
621  for (; i < n; ++i)
622  res = res * x;
623 
624  return res;
625 
626  } else if(n < 0) {
627  return T(1.0) / pow(x, -n);
628  } else {
629  return 1;
630  }
631  }
632 
633 
643  template<typename T = real>
644  TH_CONSTEXPR inline T ipow(T x, unsigned int n, T neutral_element = T(1)) {
645 
646  if(n == 0)
647  return neutral_element;
648 
649  T res = x;
650  T x_sqr = x * x;
651  unsigned int i = 1;
652 
653  // Self-multiply up to biggest power of 2
654  for (; i <= (n / 2); i *= 2)
655  res = res * res;
656 
657  // Multiply by x^2 for remaining even powers
658  for (; i < (n - 1); i += 2)
659  res = res * x_sqr;
660 
661  // Multiply for remaining powers
662  for (; i < n; ++i)
663  res = res * x;
664 
665  return res;
666  }
667 
668 
670  template<typename IntType = uint64_t>
671  TH_CONSTEXPR inline IntType fact(unsigned int n) {
672 
673  IntType res = 1;
674  for (int i = n; i > 1; --i)
675  res *= i;
676 
677  return res;
678  }
679 
680 
682  template<typename T = uint64_t>
683  TH_CONSTEXPR inline T falling_fact(T x, unsigned int n) {
684 
685  T res = 1;
686  for (unsigned int i = 0; i < n; i++)
687  res *= (x - i);
688 
689  return res;
690  }
691 
692 
694  template<typename T = uint64_t>
695  TH_CONSTEXPR inline T rising_fact(T x, unsigned int n) {
696 
697  T res = 1;
698  for (unsigned int i = 0; i < n; i++)
699  res *= (x + i);
700 
701  return res;
702  }
703 
704 
706  template<typename IntType = unsigned long long int>
707  TH_CONSTEXPR inline IntType double_fact(unsigned int n) {
708 
709  IntType res = 1;
710 
711  for (int i = n; i > 1; i -= 2)
712  res *= i;
713 
714  return res;
715  }
716 
717 
718 #ifdef THEORETICA_X86
719 
722  inline real exp_x86_norm(real x) {
723 
724  // e^x is computed as 2^(x / ln2)
725  return square(f2xm1(x / (2 * LN2)) + 1);
726  }
727 
728 #endif
729 
730 
739  inline real exp(real x) {
740 
741  // Domain reduction to [0, +inf]
742  if(x < 0)
743  return 1.0 / exp(-x);
744 
745  const real fract_x = fract(x);
746  const int floor_x = floor(x);
747 
748  // Taylor series expansion
749  // Compute e^floor(x) * e^fract(x)
750 
751  real res = 1;
752  real s_n = 1;
753 
754  for (int i = 1; i <= CORE_TAYLOR_ORDER; ++i) {
755 
756  // Recurrence formula to improve
757  // numerical stability and performance
758  s_n *= fract_x / (i * 4);
759  res += s_n;
760  }
761 
762  // The fractional part is divided by 4 to improve convergence
763  const real sqr_r = res * res;
764  return pow(E, floor_x) * sqr_r * sqr_r;
765  }
766 
767 
774  inline real expm1(real x) {
775 
776  if(abs(x) > 0.001)
777  return exp(x) - 1;
778 
779  real res = 0;
780  real s_n = 1;
781 
782  for (int i = 1; i <= 4; ++i) {
783  s_n *= x / i;
784  res += s_n;
785  }
786 
787  return res;
788  }
789 
790 
796  inline real powf(real x, real a) {
797 
798  if(a < 0)
799  return 1.0 / exp(-a * ln(abs(x)) * sgn(x));
800 
801  // x^a = e^(a * ln(x))
802  return exp(a * ln(abs(x)) * sgn(x));
803  }
804 
805 
813  inline real root(real x, int n) {
814 
815  if(((n % 2 == 0) && (x < 0)) || (n == 0)) {
816  TH_MATH_ERROR("root", n, OUT_OF_DOMAIN);
817  return nan();
818  }
819 
820  if(n < 0)
821  return 1.0 / root(x, -n);
822 
823  // Trivial cases
824  if(n == 1)
825  return x;
826 
827  if(n == 2)
828  return sqrt(x);
829 
830  if(n == 3)
831  return cbrt(x);
832 
833  if(x < 1) {
834 
835  if(x == 0)
836  return 0;
837 
838  // Approximate root between 0 and 1
839  // The root of the inverse is the inverse of the root
840  // !!! Possible precision problems with smaller numbers
841  return 1.0 / root(1.0 / x, n);
842  }
843 
844  // Approximate n-th root using Newton-Raphson.
845  // If fast exponentials and logarithms are available,
846  // use a first calculation to speed up convergence.
847 #ifdef THEORETICA_X86
848  real y = exp(ln(x) / n);
849 #else
850  real y = x;
851 #endif
852 
853  unsigned int i = 0;
854 
855  while(i < OPTIMIZATION_NEWTON_ITER) {
856 
857  const real y_pow = pow(y, n - 1);
858 
859  if((y_pow * y - x) < OPTIMIZATION_TOL)
860  break;
861 
862  y = (y * (n - 1) + x / y_pow) / (real) n;
863  i++;
864  }
865 
866  if(i >= OPTIMIZATION_NEWTON_ITER) {
867  TH_MATH_ERROR("root", i, NO_ALGO_CONVERGENCE);
868  return nan();
869  }
870 
871  return y;
872  }
873 
874 
881  inline real sin(real x) {
882 
883 #ifdef THEORETICA_X86
884 
885  #ifdef MSVC_ASM
886  __asm {
887  fld x
888  fsin
889  }
890  #else
891  asm("fsin" : "+t"(x));
892  #endif
893 
894  return x;
895 
896 #else
897 
898  // Clamp x between -2PI and 2PI
899  if(abs(x) >= TAU)
900  x -= floor(x / TAU) * TAU;
901 
902  // Domain reduction to [-PI, PI]
903  if(x > PI)
904  x = PI - x;
905  else if(x < -PI)
906  x = -PI - x;
907 
908  // Compute series with recurrence formula
909  real res = x;
910  real s = x;
911  const real sqr_x = x * x;
912 
913  for (int i = 1; i < 16; ++i) {
914  s = s * -sqr_x / (4 * i * i + 2 * i);
915  res += s;
916  }
917 
918  return res;
919 #endif
920  }
921 
922 
929  inline real cos(real x) {
930 
931 #ifdef THEORETICA_X86
932 
933  #ifdef MSVC_ASM
934  __asm {
935  fld x
936  fcos
937  }
938  #else
939  asm("fcos" : "+t"(x));
940  #endif
941 
942  return x;
943 
944 #else
945  return sin(PI2 - x);
946 #endif
947  }
948 
949 
956  inline real tan(real x) {
957 
958 #if defined(THEORETICA_X86) && !defined(MSVC_ASM)
959 
960  real s, c;
961 
962  asm ("fsincos" : "=t"(c), "=u"(s) : "0"(x));
963 
964  if(abs(c) < MACH_EPSILON) {
965  TH_MATH_ERROR("tan", c, DIV_BY_ZERO);
966  return nan();
967  }
968 
969  return s / c;
970 #else
971 
972  // Reflection
973  if(x < 0)
974  return -tan(-x);
975 
976  // Domain reduction to [0, PI]
977  x -= floor(x / PI) * PI;
978 
979  // Domain reduction to [0, PI / 4]
980  if(x > (PI / 4)) {
981  const real t = tan(x - PI / 4.0);
982  return (1 + t) / (1 - t);
983  }
984 
985  const real s = sin(x);
986  const real c = cos(x);
987 
988  if(abs(c) < MACH_EPSILON) {
989  TH_MATH_ERROR("tan", c, DIV_BY_ZERO);
990  return nan();
991  }
992 
993  return s / c;
994 #endif
995  }
996 
997 
1004  inline real cot(real x) {
1005 
1006  real s, c;
1007 
1008 #if defined(THEORETICA_X86) && !defined(MSVC_ASM)
1009 
1010  asm ("fsincos" : "=t"(c), "=u"(s) : "0"(x));
1011 #else
1012  s = sin(x);
1013  c = cos(x);
1014 #endif
1015 
1016  if(abs(s) < MACH_EPSILON) {
1017  TH_MATH_ERROR("cot", s, DIV_BY_ZERO);
1018  return nan();
1019  }
1020 
1021  return c / s;
1022  }
1023 
1024 
1032  inline real atan(real x) {
1033 
1034  // Domain reduction to [-1, 1]
1035  if(abs(x) > 1.0)
1036  return (PI2 - atan(1.0 / abs(x))) * sgn(x);
1037 
1038  const real x2 = x * x;
1039 
1040  // Interpolating Chebyshev polynomial
1041  // of degree 17
1042  return x * (0.999999981788655
1043  + x2 * (-0.3333303670928597
1044  + x2 * (0.1999187202864565
1045  + x2 * (-0.1419779780241299
1046  + x2 * (0.1061837062890163
1047  + x2 * (-0.07456854814404323
1048  + x2 * (0.04213762366862284
1049  + x2 * (-0.0157312490955519
1050  + x2 * 0.002766283502978695))))))));
1051  }
1052 
1053 
1060  inline real asin(real x) {
1061 
1062  if(abs(x) > 1) {
1063  TH_MATH_ERROR("asin", x, OUT_OF_DOMAIN);
1064  return nan();
1065  }
1066 
1067  return atan(x / sqrt(1 - x * x));
1068  }
1069 
1070 
1078  inline real acos(real x) {
1079 
1080  if(abs(x) > 1) {
1081  TH_MATH_ERROR("acos", x, OUT_OF_DOMAIN);
1082  return nan();
1083  }
1084 
1085  if(x < 0)
1086  return atan(sqrt(1 - x * x) / x) + PI;
1087  else
1088  return atan(sqrt(1 - x * x) / x);
1089  }
1090 
1091 
1099  inline real atan2(real y, real x) {
1100 
1101  if(x == 0) {
1102 
1103  if(y == 0) {
1104  TH_MATH_ERROR("atan2", y, OUT_OF_DOMAIN);
1105  return nan();
1106  }
1107 
1108  return sgn(y) * PI2;
1109  }
1110 
1111  if(x > 0) {
1112  return sgn(y) * atan(y / x);
1113  } else {
1114  return sgn(y) * atan(y / -x) + PI2;
1115  }
1116  }
1117 
1118 
1124  inline real sinh(real x) {
1125  real exp_x = exp(x);
1126  return (exp_x - 1.0 / exp_x) / 2.0;
1127  }
1128 
1129 
1135  inline real cosh(real x) {
1136  real exp_x = exp(x);
1137  return (exp_x + 1.0 / exp_x) / 2.0;
1138  }
1139 
1140 
1144  inline real tanh(real x) {
1145  real exp_x = exp(x);
1146  return (exp_x - 1.0 / exp_x) / (exp_x + 1.0 / exp_x);
1147  }
1148 
1149 
1153  inline real coth(real x) {
1154  real exp_x = exp(x);
1155  return (exp_x + 1.0 / exp_x) / (exp_x - 1.0 / exp_x);
1156  }
1157 
1158 
1160  inline real asinh(real x) {
1161  return ln(x + sqrt(square(x) + 1));
1162  }
1163 
1164 
1166  inline real acosh(real x) {
1167 
1168  if(x < 1) {
1169  TH_MATH_ERROR("acosh", x, OUT_OF_DOMAIN);
1170  return nan();
1171  }
1172 
1173  return ln(x + sqrt(square(x) - 1));
1174  }
1175 
1176 
1178  inline real atanh(real x) {
1179 
1180  if(x < -1 || x > 1) {
1181  TH_MATH_ERROR("atanh", x, OUT_OF_DOMAIN);
1182  return nan();
1183  }
1184 
1185  return 0.5 * ln((x + 1) / (1 - x));
1186  }
1187 
1188 
1193  inline real sigmoid(real x) {
1194  return 1.0 / (1.0 + 1.0 / exp(x));
1195  }
1196 
1197 
1202  inline real sinc(real x) {
1203 
1204  if(abs(x) <= MACH_EPSILON)
1205  return 1;
1206 
1207  return sin(PI * x) / (PI * x);
1208  }
1209 
1210 
1215  inline real heaviside(real x) {
1216 
1217  if(abs(x) < MACH_EPSILON)
1218  return 0.5;
1219 
1220  return x > 0 ? 1 : 0;
1221  }
1222 
1223 
1229  template<typename IntType = unsigned long long int>
1230  TH_CONSTEXPR inline IntType binomial_coeff(unsigned int n, unsigned int m) {
1231 
1232  if(n < m) {
1233  TH_MATH_ERROR("binomial_coeff", n, IMPOSSIBLE_OPERATION);
1234  return 0;
1235  }
1236 
1237  // TO-DO Check out of range
1238 
1239  IntType res = 1;
1240 
1241  for (unsigned int i = n; i > m; --i)
1242  res *= i;
1243 
1244  return res / fact(n - m);
1245  }
1246 
1247 
1254  return degrees * DEG2RAD;
1255  }
1256 
1257 
1264  return radians * RAD2DEG;
1265  }
1266 
1267 
1272  template<typename T>
1273  TH_CONSTEXPR inline T kronecker_delta(T i, T j) {
1274  return i == j ? 1 : 0;
1275  }
1276 
1277 
1279  template<typename IntType = unsigned long long int>
1280  TH_CONSTEXPR inline IntType catalan(unsigned int n) {
1281  return binomial_coeff(2 * n, n) / (n + 1);
1282  }
1283 
1284 }
1285 
1286 #endif
Mathematical constants and default algorithm parameters.
#define TH_CONSTEXPR
Enable constexpr in function declarations if C++14 is supported.
Definition: constants.h:151
Fundamental type traits.
Error handling.
#define TH_MATH_ERROR(F_NAME, VALUE, EXCEPTION)
TH_MATH_ERROR is a macro which throws exceptions or modifies errno (depending on which compiling opti...
Definition: error.h:219
Main namespace of the library which contains all functions and objects.
Definition: algebra.h:27
real heaviside(real x)
Compute the heaviside function.
Definition: real_analysis.h:1215
double real
A real number, defined as a floating point type.
Definition: constants.h:188
real acosh(real x)
Compute the inverse hyperbolic cosine.
Definition: real_analysis.h:1166
TH_CONSTEXPR IntType binomial_coeff(unsigned int n, unsigned int m)
Compute the binomial coefficient.
Definition: real_analysis.h:1230
auto min(const Vector &X)
Finds the minimum value inside a dataset.
Definition: dataset.h:351
dual2 sqrt(dual2 x)
Compute the square root of a second order dual number.
Definition: dual2_functions.h:54
dual cosh(dual x)
Compute the hyperbolic cosine of a dual number.
Definition: dual_functions.h:194
TH_CONSTEXPR T rising_fact(T x, unsigned int n)
Compute the rising factorial of n.
Definition: real_analysis.h:695
dual2 ln(dual2 x)
Compute the natural logarithm of a second order dual number.
Definition: dual2_functions.h:139
real sinc(real x)
Compute the normalized sinc function.
Definition: real_analysis.h:1202
dual sinh(dual x)
Compute the hyperbolic sine of a dual number.
Definition: dual_functions.h:186
real inf()
Return positive infinity in floating point representation.
Definition: error.h:76
constexpr real PI2
Half of Pi.
Definition: constants.h:209
TH_CONSTEXPR real radians(real degrees)
Convert degrees to radians.
Definition: real_analysis.h:1253
dual2 abs(dual2 x)
Compute the absolute value of a second order dual number.
Definition: dual2_functions.h:183
dual2 log2(dual2 x)
Compute the natural logarithm of a second order dual number.
Definition: dual2_functions.h:153
TH_CONSTEXPR IntType catalan(unsigned int n)
The n-th Catalan number.
Definition: real_analysis.h:1280
TH_CONSTEXPR T kronecker_delta(T i, T j)
Kronecker delta, equals 1 if i is equal to j, 0 otherwise.
Definition: real_analysis.h:1273
dual2 asin(dual2 x)
Compute the arcsine of a second order dual number.
Definition: dual2_functions.h:189
dual2 exp(dual2 x)
Compute the exponential of a second order dual number.
Definition: dual2_functions.h:130
constexpr real LOG210
The binary logarithm of 10.
Definition: constants.h:233
UnsignedIntType ilog2(UnsignedIntType x)
Find the integer logarithm of x.
Definition: real_analysis.h:551
dual2 log10(dual2 x)
Compute the natural logarithm of a second order dual number.
Definition: dual2_functions.h:168
TH_CONSTEXPR T ipow(T x, unsigned int n, T neutral_element=T(1))
Compute the n-th positive power of x (where n is natural)
Definition: real_analysis.h:644
real asinh(real x)
Compute the inverse hyperbolic sine.
Definition: real_analysis.h:1160
TH_CONSTEXPR IntType double_fact(unsigned int n)
Compute the double factorial of n.
Definition: real_analysis.h:707
constexpr real OPTIMIZATION_TOL
Approximation tolerance for root finding.
Definition: constants.h:278
dual2 conjugate(dual2 x)
Return the conjugate of a second order dual number.
Definition: dual2_functions.h:35
real cbrt(real x)
Compute the cubic root of x.
Definition: real_analysis.h:199
UnsignedIntType isqrt(UnsignedIntType n)
Compute the integer square root of a positive integer.
Definition: real_analysis.h:55
real atan2(real y, real x)
Compute the 2 argument arctangent.
Definition: real_analysis.h:1099
real powf(real x, real a)
Approximate x elevated to a real exponent.
Definition: real_analysis.h:796
constexpr real LN2
The natural logarithm of 2.
Definition: constants.h:239
constexpr real LOG2E
The binary logarithm of e.
Definition: constants.h:230
auto max(const Vector &X)
Finds the maximum value inside a dataset.
Definition: dataset.h:330
real sigmoid(real x)
Compute the sigmoid function.
Definition: real_analysis.h:1193
TH_CONSTEXPR real degrees(real radians)
Convert radians to degrees.
Definition: real_analysis.h:1263
real fract(real x)
Compute the fractional part of a real number.
Definition: real_analysis.h:288
TH_CONSTEXPR IntType fact(unsigned int n)
Compute the factorial of n.
Definition: real_analysis.h:671
real coth(real x)
Compute the hyperbolic cotangent.
Definition: real_analysis.h:1153
dual2 cos(dual2 x)
Compute the cosine of a second order dual number.
Definition: dual2_functions.h:82
dual2 cot(dual2 x)
Compute the cotangent of a second order dual number.
Definition: dual2_functions.h:112
constexpr real MACH_EPSILON
Machine epsilon for the real type.
Definition: constants.h:197
constexpr real E
The Euler mathematical constant (e)
Definition: constants.h:227
UnsignedIntType pad2(UnsignedIntType x)
Get the smallest power of 2 bigger than or equal to x.
Definition: real_analysis.h:578
real atanh(real x)
Compute the inverse hyperbolic tangent.
Definition: real_analysis.h:1178
constexpr int CORE_TAYLOR_ORDER
Order of Taylor series approximations.
Definition: constants.h:269
dual2 tan(dual2 x)
Compute the tangent of a second order dual number.
Definition: dual2_functions.h:94
constexpr unsigned int OPTIMIZATION_NEWTON_ITER
Maximum number of iterations for the Newton-Raphson algorithm.
Definition: constants.h:290
dual2 acos(dual2 x)
Compute the arcosine of a second order dual number.
Definition: dual2_functions.h:206
real root(real x, int n)
Compute the n-th root of x.
Definition: real_analysis.h:813
int sgn(real x)
Return the sign of x (1 if positive, -1 if negative, 0 if null)
Definition: real_analysis.h:259
dual2 sin(dual2 x)
Compute the sine of a second order dual number.
Definition: dual2_functions.h:70
UnsignedIntType icbrt(UnsignedIntType n)
Compute the integer cubic root of a positive integer.
Definition: real_analysis.h:95
TH_CONSTEXPR T falling_fact(T x, unsigned int n)
Compute the falling factorial of n.
Definition: real_analysis.h:683
complex< T > identity(complex< T > z)
Complex identity.
Definition: complex_analysis.h:19
constexpr real DEG2RAD
The scalar conversion factor from degrees to radians.
Definition: constants.h:245
real expm1(real x)
Compute the exponential of x minus 1 more accurately for really small x.
Definition: real_analysis.h:774
dual2 square(dual2 x)
Return the square of a second order dual number.
Definition: dual2_functions.h:23
constexpr real RAD2DEG
The scalar conversion factor from radians to degrees.
Definition: constants.h:248
constexpr real PI
The Pi mathematical constant.
Definition: constants.h:206
real nan()
Return a quiet NaN number in floating point representation.
Definition: error.h:54
constexpr real TAU
The Tau mathematical constant (Pi times 2)
Definition: constants.h:218
dual2 atan(dual2 x)
Compute the arctangent of a second order dual number.
Definition: dual2_functions.h:223
real clamp(real x, real a, real b)
Clamp x between a and b.
Definition: real_analysis.h:355
TH_CONSTEXPR int floor(real x)
Compute the floor of x Computes the maximum integer number that is smaller than x.
Definition: real_analysis.h:271
dual tanh(dual x)
Compute the hyperbolic tangent of a dual number.
Definition: dual_functions.h:202
dual2 pow(dual2 x, int n)
Compute the n-th power of a second order dual number.
Definition: dual2_functions.h:41
dual2 cube(dual2 x)
Return the cube of a second order dual number.
Definition: dual2_functions.h:29