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.0)
447  return -log2(1.0 / x);
448 
449  // Compute the biggest power of 2 so that x <= 2^i
450  unsigned int i = 0;
451  while(x > (uint64_t(1) << i))
452  i++;
453 
454  // Domain reduction to [1, 2]
455  x /= (uint64_t(1) << i);
456 
457  // Use the Taylor expansion of the logarithm
458  // ln(1 + x) = \sum_k^n (-1)^(k+1) x^k / k
459  real log_z = 0;
460  const real z = x - 1;
461 
462  // Exact powers of 2 don't need further computation
463  if(abs(z) > MACH_EPSILON) {
464 
465  int sign = 1;
466  real pow_z = z;
467  log_z = z;
468 
469  for (int i = 2; i <= 24; ++i) {
470 
471  pow_z *= z;
472  sign *= -1;
473 
474  log_z += sign * pow_z / i;
475  }
476  }
477 
478  return i + (log_z / LN2);
479 #endif
480  }
481 
482 
490  inline real log10(real x) {
491 
492  if(x <= 0) {
493 
494  if(x == 0) {
495  TH_MATH_ERROR("log10", x, OUT_OF_RANGE);
496  return -inf();
497  }
498 
499  TH_MATH_ERROR("log10", x, OUT_OF_DOMAIN);
500  return nan();
501  }
502 
503 #ifdef THEORETICA_X86
504 
505  // Approximate the binary logarithm of x by
506  // exploiting x86 Assembly instructions
507  return fyl2x(x, 1.0 / LOG210);
508 #else
509  return log2(x) / LOG210;
510 #endif
511  }
512 
513 
521  inline real ln(real x) {
522 
523  if(x <= 0) {
524 
525  if(x == 0) {
526  TH_MATH_ERROR("ln", x, OUT_OF_RANGE);
527  return -inf();
528  }
529 
530  TH_MATH_ERROR("ln", x, OUT_OF_DOMAIN);
531  return nan();
532  }
533 
534 #ifdef THEORETICA_X86
535 
536  // Approximate the binary logarithm of x by
537  // exploiting x86 Assembly instructions
538  return fyl2x(x, 1.0 / LOG2E);
539 #else
540  return log2(x) / LOG2E;
541 #endif
542  }
543 
544 
549  template<typename UnsignedIntType = uint64_t>
550  inline UnsignedIntType ilog2(UnsignedIntType x) {
551 
552  if(x == 0) {
553  TH_MATH_ERROR("ilog2", x, OUT_OF_RANGE);
554  return -inf();
555  }
556 
557  UnsignedIntType bit = 0;
558 
559  // Find the highest set bit
560  for (int i = (sizeof(UnsignedIntType) * 8 - 1); i > 0; --i) {
561  if(x & ((UnsignedIntType) 1 << i)) {
562  bit = i;
563  break;
564  }
565  }
566 
567  return bit;
568  }
569 
570 
576  template<typename UnsignedIntType = uint64_t>
577  inline UnsignedIntType pad2(UnsignedIntType x) {
578 
579  UnsignedIntType bit = 0;
580 
581  for (int i = (sizeof(UnsignedIntType) * 8 - 1); i > 0; --i) {
582  if(x & ((UnsignedIntType) 1 << i)) {
583 
584  // Find the highest set bit
585  bit = i;
586 
587  // Check if x is a power of 2
588  for (int j = 0; j < i; ++j)
589  if(x & ((UnsignedIntType) 1 << j))
590  return (1 << (bit + 1));
591  }
592  }
593 
594  return (1 << bit);
595  }
596 
597 
602  template<typename T = real>
603  TH_CONSTEXPR inline T pow(T x, int n) {
604 
605  if(n > 0) {
606 
607  T res = x;
608  T x_sqr = x * x;
609  int i = 1;
610 
611  // Self-multiply up to biggest power of 2
612  for (; i < (n / 2); i *= 2)
613  res = res * res;
614 
615  // Multiply by x^2 for remaining even powers
616  for (; i < (n - 1); i += 2)
617  res = res * x_sqr;
618 
619  // Multiply for remaining powers
620  for (; i < n; ++i)
621  res = res * x;
622 
623  return res;
624 
625  } else if(n < 0) {
626  return T(1.0) / pow(x, -n);
627  } else {
628  return 1;
629  }
630  }
631 
632 
642  template<typename T = real>
643  TH_CONSTEXPR inline T ipow(T x, unsigned int n, T neutral_element = T(1)) {
644 
645  if(n == 0)
646  return neutral_element;
647 
648  T res = x;
649  T x_sqr = x * x;
650  unsigned int i = 1;
651 
652  // Self-multiply up to biggest power of 2
653  for (; i <= (n / 2); i *= 2)
654  res = res * res;
655 
656  // Multiply by x^2 for remaining even powers
657  for (; i < (n - 1); i += 2)
658  res = res * x_sqr;
659 
660  // Multiply for remaining powers
661  for (; i < n; ++i)
662  res = res * x;
663 
664  return res;
665  }
666 
667 
669  template<typename IntType = uint64_t>
670  TH_CONSTEXPR inline IntType fact(unsigned int n) {
671 
672  IntType res = 1;
673  for (int i = n; i > 1; --i)
674  res *= i;
675 
676  return res;
677  }
678 
679 
681  template<typename T = uint64_t>
682  TH_CONSTEXPR inline T falling_fact(T x, unsigned int n) {
683 
684  T res = 1;
685  for (unsigned int i = 0; i < n; i++)
686  res *= (x - i);
687 
688  return res;
689  }
690 
691 
693  template<typename T = uint64_t>
694  TH_CONSTEXPR inline T rising_fact(T x, unsigned int n) {
695 
696  T res = 1;
697  for (unsigned int i = 0; i < n; i++)
698  res *= (x + i);
699 
700  return res;
701  }
702 
703 
705  template<typename IntType = unsigned long long int>
706  TH_CONSTEXPR inline IntType double_fact(unsigned int n) {
707 
708  IntType res = 1;
709 
710  for (int i = n; i > 1; i -= 2)
711  res *= i;
712 
713  return res;
714  }
715 
716 
717 #ifdef THEORETICA_X86
718 
721  inline real exp_x86_norm(real x) {
722 
723  // e^x is computed as 2^(x / ln2)
724  return square(f2xm1(x / (2 * LN2)) + 1);
725  }
726 
727 #endif
728 
729 
738  inline real exp(real x) {
739 
740  // Domain reduction to [0, +inf]
741  if(x < 0)
742  return 1.0 / exp(-x);
743 
744  const real fract_x = fract(x);
745  const int floor_x = floor(x);
746 
747  // Taylor series expansion
748  // Compute e^floor(x) * e^fract(x)
749 
750  real res = 1;
751  real s_n = 1;
752 
753  for (int i = 1; i <= CORE_TAYLOR_ORDER; ++i) {
754 
755  // Recurrence formula to improve
756  // numerical stability and performance
757  s_n *= fract_x / (i * 4);
758  res += s_n;
759  }
760 
761  // The fractional part is divided by 4 to improve convergence
762  const real sqr_r = res * res;
763  return pow(E, floor_x) * sqr_r * sqr_r;
764  }
765 
766 
773  inline real expm1(real x) {
774 
775  if(abs(x) > 0.001)
776  return exp(x) - 1;
777 
778  real res = 0;
779  real s_n = 1;
780 
781  for (int i = 1; i <= 4; ++i) {
782  s_n *= x / i;
783  res += s_n;
784  }
785 
786  return res;
787  }
788 
789 
795  inline real powf(real x, real a) {
796 
797  if(a < 0)
798  return 1.0 / exp(-a * ln(abs(x)) * sgn(x));
799 
800  // x^a = e^(a * ln(x))
801  return exp(a * ln(abs(x)) * sgn(x));
802  }
803 
804 
812  inline real root(real x, int n) {
813 
814  if(((n % 2 == 0) && (x < 0)) || (n == 0)) {
815  TH_MATH_ERROR("root", n, OUT_OF_DOMAIN);
816  return nan();
817  }
818 
819  if(n < 0)
820  return 1.0 / root(x, -n);
821 
822  // Trivial cases
823  if(n == 1)
824  return x;
825 
826  if(n == 2)
827  return sqrt(x);
828 
829  if(n == 3)
830  return cbrt(x);
831 
832  if(x < 1) {
833 
834  if(x == 0)
835  return 0;
836 
837  // Approximate root between 0 and 1
838  // The root of the inverse is the inverse of the root
839  // !!! Possible precision problems with smaller numbers
840  return 1.0 / root(1.0 / x, n);
841  }
842 
843  // Approximate n-th root using Newton-Raphson.
844  // If fast exponentials and logarithms are available,
845  // use a first calculation to speed up convergence.
846 #ifdef THEORETICA_X86
847  real y = exp(ln(x) / n);
848 #else
849  real y = x;
850 #endif
851 
852  unsigned int i = 0;
853 
854  while(i < OPTIMIZATION_NEWTON_ITER) {
855 
856  const real y_pow = pow(y, n - 1);
857 
858  if((y_pow * y - x) < OPTIMIZATION_TOL)
859  break;
860 
861  y = (y * (n - 1) + x / y_pow) / (real) n;
862  i++;
863  }
864 
865  if(i >= OPTIMIZATION_NEWTON_ITER) {
866  TH_MATH_ERROR("root", i, NO_ALGO_CONVERGENCE);
867  return nan();
868  }
869 
870  return y;
871  }
872 
873 
880  inline real sin(real x) {
881 
882 #ifdef THEORETICA_X86
883 
884  #ifdef MSVC_ASM
885  __asm {
886  fld x
887  fsin
888  }
889  #else
890  asm("fsin" : "+t"(x));
891  #endif
892 
893  return x;
894 
895 #else
896 
897  // Clamp x between -2PI and 2PI
898  if(abs(x) >= TAU)
899  x -= floor(x / TAU) * TAU;
900 
901  // Domain reduction to [-PI, PI]
902  if(x > PI)
903  x = PI - x;
904  else if(x < -PI)
905  x = -PI - x;
906 
907  // Compute series with recurrence formula
908  real res = x;
909  real s = x;
910  const real sqr_x = x * x;
911 
912  for (int i = 1; i < 16; ++i) {
913  s = s * -sqr_x / (4 * i * i + 2 * i);
914  res += s;
915  }
916 
917  return res;
918 #endif
919  }
920 
921 
928  inline real cos(real x) {
929 
930 #ifdef THEORETICA_X86
931 
932  #ifdef MSVC_ASM
933  __asm {
934  fld x
935  fcos
936  }
937  #else
938  asm("fcos" : "+t"(x));
939  #endif
940 
941  return x;
942 
943 #else
944  return sin(PI2 - x);
945 #endif
946  }
947 
948 
955  inline real tan(real x) {
956 
957 #if defined(THEORETICA_X86) && !defined(MSVC_ASM)
958 
959  real s, c;
960 
961  asm ("fsincos" : "=t"(c), "=u"(s) : "0"(x));
962 
963  if(abs(c) < MACH_EPSILON) {
964  TH_MATH_ERROR("tan", c, DIV_BY_ZERO);
965  return nan();
966  }
967 
968  return s / c;
969 #else
970 
971  // Reflection
972  if(x < 0)
973  return -tan(-x);
974 
975  // Domain reduction to [0, PI]
976  x -= floor(x / PI) * PI;
977 
978  // Domain reduction to [0, PI / 4]
979  if(x > (PI / 4)) {
980  const real t = tan(x - PI / 4.0);
981  return (1 + t) / (1 - t);
982  }
983 
984  const real s = sin(x);
985  const real c = cos(x);
986 
987  if(abs(c) < MACH_EPSILON) {
988  TH_MATH_ERROR("tan", c, DIV_BY_ZERO);
989  return nan();
990  }
991 
992  return s / c;
993 #endif
994  }
995 
996 
1003  inline real cot(real x) {
1004 
1005  real s, c;
1006 
1007 #if defined(THEORETICA_X86) && !defined(MSVC_ASM)
1008 
1009  asm ("fsincos" : "=t"(c), "=u"(s) : "0"(x));
1010 #else
1011  s = sin(x);
1012  c = cos(x);
1013 #endif
1014 
1015  if(abs(s) < MACH_EPSILON) {
1016  TH_MATH_ERROR("cot", s, DIV_BY_ZERO);
1017  return nan();
1018  }
1019 
1020  return c / s;
1021  }
1022 
1023 
1031  inline real atan(real x) {
1032 
1033  // Domain reduction to [-1, 1]
1034  if(abs(x) > 1.0)
1035  return (PI2 - atan(1.0 / abs(x))) * sgn(x);
1036 
1037  const real x2 = x * x;
1038 
1039  // Interpolating Chebyshev polynomial
1040  // of degree 17
1041  return x * (0.999999981788655
1042  + x2 * (-0.3333303670928597
1043  + x2 * (0.1999187202864565
1044  + x2 * (-0.1419779780241299
1045  + x2 * (0.1061837062890163
1046  + x2 * (-0.07456854814404323
1047  + x2 * (0.04213762366862284
1048  + x2 * (-0.0157312490955519
1049  + x2 * 0.002766283502978695))))))));
1050  }
1051 
1052 
1059  inline real asin(real x) {
1060 
1061  if(abs(x) > 1) {
1062  TH_MATH_ERROR("asin", x, OUT_OF_DOMAIN);
1063  return nan();
1064  }
1065 
1066  return atan(x / sqrt(1 - x * x));
1067  }
1068 
1069 
1077  inline real acos(real x) {
1078 
1079  if(abs(x) > 1) {
1080  TH_MATH_ERROR("acos", x, OUT_OF_DOMAIN);
1081  return nan();
1082  }
1083 
1084  if(x < 0)
1085  return atan(sqrt(1 - x * x) / x) + PI;
1086  else
1087  return atan(sqrt(1 - x * x) / x);
1088  }
1089 
1090 
1098  inline real atan2(real y, real x) {
1099 
1100  if(x == 0) {
1101 
1102  if(y == 0) {
1103  TH_MATH_ERROR("atan2", y, OUT_OF_DOMAIN);
1104  return nan();
1105  }
1106 
1107  return sgn(y) * PI2;
1108  }
1109 
1110  if(x > 0) {
1111  return sgn(y) * atan(y / x);
1112  } else {
1113  return sgn(y) * atan(y / -x) + PI2;
1114  }
1115  }
1116 
1117 
1123  inline real sinh(real x) {
1124  real exp_x = exp(x);
1125  return (exp_x - 1.0 / exp_x) / 2.0;
1126  }
1127 
1128 
1134  inline real cosh(real x) {
1135  real exp_x = exp(x);
1136  return (exp_x + 1.0 / exp_x) / 2.0;
1137  }
1138 
1139 
1143  inline real tanh(real x) {
1144  real exp_x = exp(x);
1145  return (exp_x - 1.0 / exp_x) / (exp_x + 1.0 / exp_x);
1146  }
1147 
1148 
1152  inline real coth(real x) {
1153  real exp_x = exp(x);
1154  return (exp_x + 1.0 / exp_x) / (exp_x - 1.0 / exp_x);
1155  }
1156 
1157 
1159  inline real asinh(real x) {
1160  return ln(x + sqrt(square(x) + 1));
1161  }
1162 
1163 
1165  inline real acosh(real x) {
1166 
1167  if(x < 1) {
1168  TH_MATH_ERROR("acosh", x, OUT_OF_DOMAIN);
1169  return nan();
1170  }
1171 
1172  return ln(x + sqrt(square(x) - 1));
1173  }
1174 
1175 
1177  inline real atanh(real x) {
1178 
1179  if(x < -1 || x > 1) {
1180  TH_MATH_ERROR("atanh", x, OUT_OF_DOMAIN);
1181  return nan();
1182  }
1183 
1184  return 0.5 * ln((x + 1) / (1 - x));
1185  }
1186 
1187 
1192  inline real sigmoid(real x) {
1193  return 1.0 / (1.0 + 1.0 / exp(x));
1194  }
1195 
1196 
1201  inline real sinc(real x) {
1202 
1203  if(abs(x) <= MACH_EPSILON)
1204  return 1;
1205 
1206  return sin(PI * x) / (PI * x);
1207  }
1208 
1209 
1214  inline real heaviside(real x) {
1215 
1216  if(abs(x) < MACH_EPSILON)
1217  return 0.5;
1218 
1219  return x > 0 ? 1 : 0;
1220  }
1221 
1222 
1228  template<typename IntType = unsigned long long int>
1229  TH_CONSTEXPR inline IntType binomial_coeff(unsigned int n, unsigned int m) {
1230 
1231  if(n < m) {
1232  TH_MATH_ERROR("binomial_coeff", n, IMPOSSIBLE_OPERATION);
1233  return 0;
1234  }
1235 
1236  // TO-DO Check out of range
1237 
1238  IntType res = 1;
1239 
1240  for (unsigned int i = n; i > m; --i)
1241  res *= i;
1242 
1243  return res / fact(n - m);
1244  }
1245 
1246 
1253  return degrees * DEG2RAD;
1254  }
1255 
1256 
1263  return radians * RAD2DEG;
1264  }
1265 
1266 
1271  template<typename T>
1272  TH_CONSTEXPR inline T kronecker_delta(T i, T j) {
1273  return i == j ? 1 : 0;
1274  }
1275 
1276 
1278  template<typename IntType = unsigned long long int>
1279  TH_CONSTEXPR inline IntType catalan(unsigned int n) {
1280  return binomial_coeff(2 * n, n) / (n + 1);
1281  }
1282 
1283 }
1284 
1285 #endif
Mathematical constants and default algorithm parameters.
#define TH_CONSTEXPR
Enable constexpr in function declarations if C++14 is supported.
Definition: constants.h:161
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:1214
double real
A real number, defined as a floating point type.
Definition: constants.h:198
real acosh(real x)
Compute the inverse hyperbolic cosine.
Definition: real_analysis.h:1165
TH_CONSTEXPR IntType binomial_coeff(unsigned int n, unsigned int m)
Compute the binomial coefficient.
Definition: real_analysis.h:1229
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:205
TH_CONSTEXPR T rising_fact(T x, unsigned int n)
Compute the rising factorial of n.
Definition: real_analysis.h:694
dual2 ln(dual2 x)
Compute the natural logarithm of a second order dual number.
Definition: dual2_functions.h:151
real sinc(real x)
Compute the normalized sinc function.
Definition: real_analysis.h:1201
dual sinh(dual x)
Compute the hyperbolic sine of a dual number.
Definition: dual_functions.h:194
real inf()
Return positive infinity in floating point representation.
Definition: error.h:76
constexpr real PI2
Half of Pi.
Definition: constants.h:219
TH_CONSTEXPR real radians(real degrees)
Convert degrees to radians.
Definition: real_analysis.h:1252
dual2 abs(dual2 x)
Compute the absolute value of a second order dual number.
Definition: dual2_functions.h:198
dual2 log2(dual2 x)
Compute the natural logarithm of a second order dual number.
Definition: dual2_functions.h:166
TH_CONSTEXPR IntType catalan(unsigned int n)
The n-th Catalan number.
Definition: real_analysis.h:1279
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:1272
dual2 asin(dual2 x)
Compute the arcsine of a second order dual number.
Definition: dual2_functions.h:204
dual2 exp(dual2 x)
Compute the exponential of a second order dual number.
Definition: dual2_functions.h:138
constexpr real LOG210
The binary logarithm of 10.
Definition: constants.h:243
UnsignedIntType ilog2(UnsignedIntType x)
Find the integer logarithm of x.
Definition: real_analysis.h:550
dual2 log10(dual2 x)
Compute the natural logarithm of a second order dual number.
Definition: dual2_functions.h:182
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:643
real asinh(real x)
Compute the inverse hyperbolic sine.
Definition: real_analysis.h:1159
TH_CONSTEXPR IntType double_fact(unsigned int n)
Compute the double factorial of n.
Definition: real_analysis.h:706
constexpr real OPTIMIZATION_TOL
Approximation tolerance for root finding.
Definition: constants.h:288
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:1098
real powf(real x, real a)
Approximate x elevated to a real exponent.
Definition: real_analysis.h:795
constexpr real LN2
The natural logarithm of 2.
Definition: constants.h:249
constexpr real LOG2E
The binary logarithm of e.
Definition: constants.h:240
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:1192
TH_CONSTEXPR real degrees(real radians)
Convert radians to degrees.
Definition: real_analysis.h:1262
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:670
real coth(real x)
Compute the hyperbolic cotangent.
Definition: real_analysis.h:1152
dual2 cos(dual2 x)
Compute the cosine of a second order dual number.
Definition: dual2_functions.h:86
dual2 cot(dual2 x)
Compute the cotangent of a second order dual number.
Definition: dual2_functions.h:119
constexpr real MACH_EPSILON
Machine epsilon for the real type.
Definition: constants.h:207
constexpr real E
The Euler mathematical constant (e)
Definition: constants.h:237
UnsignedIntType pad2(UnsignedIntType x)
Get the smallest power of 2 bigger than or equal to x.
Definition: real_analysis.h:577
real atanh(real x)
Compute the inverse hyperbolic tangent.
Definition: real_analysis.h:1177
constexpr int CORE_TAYLOR_ORDER
Order of Taylor series approximations.
Definition: constants.h:279
dual2 tan(dual2 x)
Compute the tangent of a second order dual number.
Definition: dual2_functions.h:100
constexpr unsigned int OPTIMIZATION_NEWTON_ITER
Maximum number of iterations for the Newton-Raphson algorithm.
Definition: constants.h:300
dual2 acos(dual2 x)
Compute the arcosine of a second order dual number.
Definition: dual2_functions.h:223
real root(real x, int n)
Compute the n-th root of x.
Definition: real_analysis.h:812
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:72
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:682
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:255
real expm1(real x)
Compute the exponential of x minus 1 more accurately for really small x.
Definition: real_analysis.h:773
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:258
constexpr real PI
The Pi mathematical constant.
Definition: constants.h:216
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:228
dual2 atan(dual2 x)
Compute the arctangent of a second order dual number.
Definition: dual2_functions.h:242
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:216
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