/usr/include/scythestat/smath.h is in libscythestat-dev 1.0.2-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 | /*
* Scythe Statistical Library Copyright (C) 2000-2002 Andrew D. Martin
* and Kevin M. Quinn; 2002-present Andrew D. Martin, Kevin M. Quinn,
* and Daniel Pemstein. All Rights Reserved.
*
* This program is free software; you can redistribute it and/or
* modify under the terms of the GNU General Public License as
* published by Free Software Foundation; either version 2 of the
* License, or (at your option) any later version. See the text files
* COPYING and LICENSE, distributed with this source code, for further
* information.
* --------------------------------------------------------------------
* scythestat/smath.h
*
*/
/*!
* \file smath.h
* \brief Definitions for functions that perform common mathematical
* operations on every element of a Matrix.
*
* \note As is the case throughout the library, we provide both
* general and default template definitions of the Matrix-returning
* functions in this file, explicitly providing documentation for only
* the general template versions. As is also often the case, Doxygen
* does not always correctly add the default template definition to
* the function list below; there is always a default template
* definition available for every function.
*
*/
#ifndef SCYTHE_MATH_H
#define SCYTHE_MATH_H
#ifdef SCYTHE_COMPILE_DIRECT
#include "matrix.h"
#include "algorithm.h"
#include "error.h"
#else
#include "scythestat/matrix.h"
#include "scythestat/algorithm.h"
#include "scythestat/error.h"
#endif
#include <cmath>
#include <numeric>
#include <set>
namespace scythe {
namespace {
typedef unsigned int uint;
}
/* Almost every function in this file follows one of the two patterns
* described by these macros. The first macro handles single-argument
* functions. The second handles two-matrix-argument functions (or
* scalar-matrix, matrix-scalar. The second macro also permits
* cross-type operations (these are limited only by the capabilities
* of the underlying functions).
*/
#define SCYTHE_MATH_OP(NAME, OP) \
template <matrix_order RO, matrix_style RS, typename T, \
matrix_order PO, matrix_style PS> \
Matrix<T,RO,RS> \
NAME (const Matrix<T,PO,PS>& A) \
{ \
Matrix<T,RO,RS> res(A.rows(), A.cols(), false); \
std::transform(A.begin_f(), A.end_f(), res.begin_f(), OP); \
return res; \
} \
\
template <typename T, matrix_order O, matrix_style S> \
Matrix<T,O,Concrete> \
NAME (const Matrix<T,O,S>& A) \
{ \
return NAME<O,Concrete>(A); \
}
#define SCYTHE_MATH_OP_2ARG(NAME, OP) \
template <matrix_order RO, matrix_style RS, typename T, \
matrix_order PO1, matrix_style PS1, \
matrix_order PO2, matrix_style PS2, typename S> \
Matrix<T,RO,RS> \
NAME (const Matrix<T,PO1,PS1>& A, const Matrix<S,PO2,PS2>& B) \
{ \
SCYTHE_CHECK_10 (A.size() != 1 && B.size() != 1 && \
A.size() != B.size(), scythe_conformation_error, \
"Matrices with dimensions (" << A.rows() \
<< ", " << A.cols() \
<< ") and (" << B.rows() << ", " << B.cols() \
<< ") are not conformable"); \
\
Matrix<T,RO,RS> res; \
\
if (A.size() == 1) { \
res.resize2Match(B); \
std::transform(B.template begin_f<RO>(), B.template end_f<RO>(),\
res.begin_f(), std::bind1st(OP, A(0))); \
} else if (B.size() == 1) { \
res.resize2Match(A); \
std::transform(A.template begin_f<RO>(), A.template end_f<RO>(),\
res.begin_f(), std::bind2nd(OP, B(0))); \
} else { \
res.resize2Match(A); \
std::transform(A.template begin_f<RO>(), A.template end_f<RO>(),\
B.template begin_f<RO>(), res.begin_f(), OP); \
} \
\
return res; \
} \
\
template <typename T, matrix_order PO1, matrix_style PS1, \
matrix_order PO2, matrix_style PS2, \
typename S> \
Matrix<T,PO1,Concrete> \
NAME (const Matrix<T,PO1,PS1>& A, const Matrix<S,PO2,PS2>& B) \
{ \
return NAME<PO1,Concrete>(A, B); \
} \
\
template<matrix_order RO, matrix_style RS, typename T, \
matrix_order PO, matrix_style PS, typename S> \
Matrix<T,RO,RS> \
NAME (const Matrix<T,PO,PS>& A, S b) \
{ \
return NAME<RO,RS>(A, Matrix<S,RO,Concrete>(b)); \
} \
\
template <typename T, typename S, matrix_order PO, matrix_style PS> \
Matrix<T,PO,Concrete> \
NAME (const Matrix<T,PO,PS>& A, S b) \
{ \
return NAME<PO,Concrete>(A, Matrix<S,PO,Concrete>(b)); \
} \
\
template<matrix_order RO, matrix_style RS, typename T, \
matrix_order PO, matrix_style PS, typename S> \
Matrix<T,RO,RS> \
NAME (T a, const Matrix<S,PO,PS>& B) \
{ \
return NAME<RO,RS>(Matrix<S, RO,Concrete>(a), B); \
} \
\
template <typename T, typename S, matrix_order PO, matrix_style PS> \
Matrix<T,PO,Concrete> \
NAME (T a, const Matrix<S,PO,PS>& B) \
{ \
return NAME<PO,Concrete>(Matrix<S,PO,Concrete>(a), B); \
}
/* calc the inverse cosine of each element of a Matrix */
/*!
* \brief Calculate the inverse cosine of each element of a Matrix
*
* This function calculates the inverse cosine of each element in a Matrix
*
* \param A The matrix whose inverse cosines are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cos()
* \see cosh()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(acos, ::acos)
/* calc the inverse hyperbolic cosine of each element of a Matrix */
/*!
* \brief Calculate the inverse hyperbolic cosine of each element of a Matrix
*
* This function calculates the inverse hyperbolic cosine of each element
* in a Matrix
*
* \param A The matrix whose inverse hyperbolic cosines are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cos()
* \see cosh()
* \see acos()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(acosh, ::acosh)
/* calc the inverse sine of each element of a Matrix */
/*!
* \brief Calculate the inverse sine of each element of a Matrix
*
* This function calculates the inverse sine of each element
* in a Matrix
*
* \param A The matrix whose inverse sines are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(asin, ::asin)
/* calc the inverse hyperbolic sine of each element of a Matrix */
/*!
* \brief Calculate the inverse hyperbolic sine of each element of a Matrix
*
* This function calculates the inverse hyperbolic sine of each element
* in a Matrix
*
* \param A The matrix whose inverse hyperbolic sines are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(asinh, ::asinh)
/* calc the inverse tangent of each element of a Matrix */
/*!
* \brief Calculate the inverse tangent of each element of a Matrix
*
* This function calculates the inverse tangent of each element
* in a Matrix
*
* \param A The matrix whose inverse tangents are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see asin()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(atan, ::atan)
/* calc the inverse hyperbolic tangent of each element of a Matrix */
/*!
* \brief Calculate the inverse hyperbolic tangent of each element of a Matrix
*
* This function calculates the inverse hyperbolic tangent of each element
* in a Matrix
*
* \param A The matrix whose inverse hyperbolic tangents are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atan2()
*/
SCYTHE_MATH_OP(atanh, ::atanh)
/* calc the angle whose tangent is y/x */
/*!
* \brief Calculate the angle whose tangent is y/x
*
* This function calculates the angle whose tangent is y/x, given two
* matrices A and B (where y is the ith element of A, and x is the jth element
* of matrix B).
*
* \param A The matrix of y values
* \param B The matrix of x values
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
*/
SCYTHE_MATH_OP_2ARG(atan2, std::ptr_fun(::atan2))
/* calc the cube root of each element of a Matrix */
/*!
* \brief Calculate the cube root of each element of a Matrix
*
* This function calculates the cube root of each element
* in a Matrix
*
* \param A The matrix whose cube roots are of interest.
*
* \see sqrt()
*/
SCYTHE_MATH_OP(cbrt, ::cbrt)
/* calc the ceil of each element of a Matrix */
/*!
* \brief Calculate the ceiling of each element of a Matrix
*
* This function calculates the ceiling of each element
* in a Matrix
*
* \param A The matrix whose ceilings are of interest.
*
* \see floor()
*/
SCYTHE_MATH_OP(ceil, ::ceil)
/* create a matrix containing the absval of the first input and the
* sign of the second
*/
/*!
* \brief Create a matrix containing the absolute value of the first input
* and the sign of the second input
*
* This function creates a matrix containing the absolute value of the first
* input, a matrix called A, and the sign of the second input, matrix B.
*
* \param A The matrix whose absolute values will comprise the resultant matrix.
* \param B The matrix whose signs will comprise the resultant matrix
*/
SCYTHE_MATH_OP_2ARG(copysign, std::ptr_fun(::copysign))
/* calc the cosine of each element of a Matrix */
/*!
* \brief Calculate the cosine of each element of a Matrix
*
* This function calculates the cosine of each element in a Matrix
*
* \param A The matrix whose cosines are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(cos, ::cos)
/* calc the hyperbolic cosine of each element of a Matrix */
/*!
* \brief Calculate the hyperbolic cosine of each element of a Matrix
*
* This function calculates the hyperbolic cosine of each element in a Matrix
*
* \param A The matrix whose hyperbolic cosines are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see sinh()
* \see cos()
* \see acos()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(cosh, ::cosh)
/* calc the error function of each element of a Matrix */
/*!
* \brief Calculate the error function of each element of a Matrix
*
* This function calculates the error function of each element in a Matrix
*
* \param A The matrix whose error functions are of interest.
*
* \see erfc()
*/
SCYTHE_MATH_OP(erf, ::erf)
/* calc the complementary error function of each element of a Matrix */
/*!
* \brief Calculate the complementary error function of each element of a Matrix
*
* This function calculates the complemenatry error function of each
* element in a Matrix
*
* \param A The matrix whose complementary error functions are of interest.
*
* \see erf()
*/
SCYTHE_MATH_OP(erfc, ::erfc)
/* calc the vaue e^x of each element of a Matrix */
/*!
* \brief Calculate the value e^x for each element of a Matrix
*
* This function calculates the value e^x for each element of a matrix, where
* x is the ith element of the matrix A
*
* \param A The matrix whose elements are to be exponentiated.
*
* \see expm1()
*/
SCYTHE_MATH_OP(exp, ::exp)
/* calc the exponent - 1 of each element of a Matrix */
/*!
* \brief Calculate the value e^(x-1) for each element of a Matrix
*
* This function calculates the value e^(x-1) for each element of a matrix, where
* x is the ith element of the matrix A
*
* \param A The matrix whose elements are to be exponentiated.
*
* \see exp()
*/
SCYTHE_MATH_OP(expm1, ::expm1)
/* calc the absval of each element of a Matrix */
/*!
* \brief Calculate the absolute value of each element of a Matrix
*
* This function calculates the absolute value of each element in a Matrix
*
* \param A The matrix whose absolute values are to be taken.
*/
SCYTHE_MATH_OP(fabs, ::fabs)
/* calc the floor of each element of a Matrix */
/*!
* \brief Calculate the floor of each element of a Matrix
*
* This function calculates the floor of each element
* in a Matrix
*
* \param A The matrix whose floors are of interest.
*
* \see ceil()
*/
SCYTHE_MATH_OP(floor, ::floor)
/* calc the remainder of the division of each matrix element */
/*!
* \brief Calculate the remainder of the division of each matrix element
*
* This function calculates the remainder when the elements of Matrix A are
* divided by the elements of Matrix B.
*
* \param A The matrix to serve as dividend
* \param B the matrix to serve as divisor
*/
SCYTHE_MATH_OP_2ARG(fmod, std::ptr_fun(::fmod))
/* calc the fractional val of input and return exponents in int
* matrix reference
*/
/*!
*/
template <matrix_order RO, matrix_style RS, typename T,
matrix_order PO1, matrix_style PS1,
matrix_order PO2, matrix_style PS2>
Matrix<T,RO,RS>
frexp (const Matrix<T,PO1,PS1>& A, Matrix<int,PO2,PS2>& ex)
{
SCYTHE_CHECK_10(A.size() != ex.size(), scythe_conformation_error,
"The input matrix sizes do not match");
Matrix<T,PO1,Concrete> res(A.rows(), A.cols());
typename Matrix<T,PO1,PS1>::const_forward_iterator it;
typename Matrix<T,PO1,Concrete>::forward_iterator rit
= res.begin_f();
typename Matrix<int,PO2,PS2>::const_forward_iterator it2
= ex.begin_f();
for (it = A.begin_f(); it != A.end_f(); ++it) {
*rit = ::frexp(*it, &(*it2));
++it2; ++rit;
}
return res;
}
template <typename T, matrix_order PO1, matrix_style PS1,
matrix_order PO2, matrix_style PS2>
Matrix<T,PO1,Concrete>
frexp (Matrix<T,PO1,PS1>& A, Matrix<int,PO2,PS2>& ex)
{
return frexp<PO1,Concrete>(A,ex);
}
/* calc the euclidean distance between the two inputs */
/*!
* \brief Calculate the euclidean distance between two inputs
*
* This function calculates the euclidean distance between the elements of Matrix
* A and the elements of Matrix B.
*
* \param A Input matrix
* \param B Input matrix
*/
SCYTHE_MATH_OP_2ARG(hypot, std::ptr_fun(::hypot))
/* return (int) logb */
SCYTHE_MATH_OP(ilogb, ::ilogb)
/* compute the bessel func of the first kind of the order 0 */
/*!
* \brief Compute the Bessel function of the first kind of the order 0
*
* This function computes the Bessel function of the first kind of order 0
* for each element in the input matrix, A.
*
* \param A Matrix for which the Bessel function is of interest
*
* \see j1()
* \see jn()
* \see y0()
* \see y1()
* \see yn()
*/
SCYTHE_MATH_OP(j0, ::j0)
/* compute the bessel func of the first kind of the order 1 */
/*!
* \brief Compute the Bessel function of the first kind of the order 1
*
* This function computes the Bessel function of the first kind of order 1
* for each element in the input matrix, A.
*
* \param A Matrix for which the Bessel function is of interest
*
* \see j0()
* \see jn()
* \see y0()
* \see y1()
* \see yn()
*/
SCYTHE_MATH_OP(j1, ::j1)
/* compute the bessel func of the first kind of the order n
* TODO: This definition causes the compiler to issue some warnings.
* Fix
*/
/*!
* \brief Compute the Bessel function of the first kind of the order n
*
* This function computes the Bessel function of the first kind of order n
* for each element in the input matrix, A.
*
* \param n Order of the Bessel function
* \param A Matrix for which the Bessel function is of interest
*
* \see j0()
* \see j1()
* \see y0()
* \see y1()
* \see yn()
*/
SCYTHE_MATH_OP_2ARG(jn, std::ptr_fun(::jn))
/* calc x * 2 ^ex */
/*!
* \brief Compute x * 2^ex
*
* This function computes the value of x * 2^ex, where x is the ith element of
* the input matrix A, and ex is the desired value of the exponent.
*
* \param A Matrix whose elements are to be multiplied
* \param ex Matrix of powers to which 2 will be raised.
*/
SCYTHE_MATH_OP_2ARG(ldexp, std::ptr_fun(::ldexp))
/* compute the natural log of the absval of gamma function */
/*!
* \brief Compute the natural log of the absolute value of the gamma function
*
* This function computes the absolute value of the Gamma Function, evaluated at
* each element of the input matrix A.
*
* \param A Matrix whose elements will serve as inputs for the Gamma Function
*
* \see log()
*/
SCYTHE_MATH_OP(lgamma, ::lgamma)
/* calc the natural log of each element of a Matrix */
/*!
* \brief Compute the natural log of each element of a Matrix
*
* This function computes the natural log of each element in a matrix, A.
*
* \param A Matrix whose natural logs are of interest
*
* \see log10()
* \see log1p()
* \see logb()
*/
SCYTHE_MATH_OP(log, ::log)
/* calc the base-10 log of each element of a Matrix */
/*!
* \brief Compute the log base 10 of each element of a Matrix
*
* This function computes the log base 10 of each element in a matrix, A.
*
* \param A Matrix whose logs are of interest
*
* \see log()
* \see log1p()
* \see logb()
*/
SCYTHE_MATH_OP(log10, ::log10)
/* calc the natural log of 1 + each element of a Matrix */
/*!
* \brief Compute the natural log of 1 + each element of a Matrix
*
* This function computes the natural log of 1 + each element of a Matrix.
*
* \param A Matrix whose logs are of interest
*
* \see log()
* \see log10()
* \see logb()
*/
SCYTHE_MATH_OP(log1p, ::log1p)
/* calc the logb of each element of a Matrix */
/*!
* \brief Compute the logb each element of a Matrix
*
* This function computes the log base b of each element of a Matrix.
*
* \param A Matrix whose logs are of interest
*
* \see log()
* \see log10()
* \see log1p()
*/
SCYTHE_MATH_OP(logb, ::logb)
/* x = frac + i, return matrix of frac and place i in 2nd matrix
*/
template <matrix_order RO, matrix_style RS, typename T,
matrix_order PO1, matrix_style PS1,
matrix_order PO2, matrix_style PS2>
Matrix<T,RO,RS>
modf (const Matrix<T,PO1,PS1>& A, Matrix<double,PO2,PS2>& ipart)
{
SCYTHE_CHECK_10(A.size() != ipart.size(), scythe_conformation_error,
"The input matrix sizes do not match");
Matrix<T,PO1,Concrete> res(A.rows(), A.cols());
typename Matrix<T,PO1,PS1>::const_forward_iterator it;
typename Matrix<T,PO1,Concrete>::forward_iterator rit
= res.begin_f();
typename Matrix<double,PO2,PS2>::const_forward_iterator it2
= ipart.begin_f();
for (it = A.begin_f(); it != A.end_f(); ++it) {
*rit = ::modf(*it, &(*it2));
++it2; ++rit;
}
return res;
}
template <typename T, matrix_order PO1, matrix_style PS1,
matrix_order PO2, matrix_style PS2>
Matrix<T,PO1,Concrete>
modf (Matrix<T,PO1,PS1>& A, Matrix<double,PO2,PS2>& ipart)
{
return modf<PO1,Concrete>(A,ipart);
}
/* calc x^ex of each element of a Matrix */
/*!
* \brief Compute x^ex for each element of a matrix
*
* This function computes x^ex, where x is the ith element of the matrix A,
* and ex is the desired exponent.
*
* \param A Matrix to be exponentiated
* \param ex Desired exponent
*/
SCYTHE_MATH_OP_2ARG(pow, std::ptr_fun(::pow))
/* calc rem == x - n * y */
SCYTHE_MATH_OP_2ARG(remainder, std::ptr_fun(::remainder))
/* return x rounded to nearest int */
/*!
* \brief Return x rounded to the nearest integer
*
* This function returns x, where x is the ith element of the Matrix A,
* rounded to the nearest integer.
*
* \param A Matrix whose elements are to be rounded
*/
SCYTHE_MATH_OP(rint, ::rint)
/* returns x * FLT_RADIX^ex */
SCYTHE_MATH_OP_2ARG(scalbn, std::ptr_fun(::scalbn))
/* calc the sine of x */
/*!
* \brief Calculate the sine of each element of a Matrix
*
* This function calculates the sine of each element in a Matrix
*
* \param A The matrix whose sines are of interest.
*
* \see tan()
* \see tanh()
* \see sinh()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(sin, ::sin)
/* calc the hyperbolic sine of x */
/*!
* \brief Calculate the hyperbolic sine of each element of a Matrix
*
* This function calculates the hyperbolic sine of each element in a Matrix
*
* \param A The matrix whose hyperbolic sines are of interest.
*
* \see tan()
* \see tanh()
* \see sin()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(sinh, ::sinh)
/* calc the sqrt of x */
/*!
* \brief Calculate the square root of each element in a matrix
*
* This function calculates the square root of each element in a Matrix
*
* \param A The matrix whose roots are of interest.
*
* \see cbrt()
*/
SCYTHE_MATH_OP(sqrt, ::sqrt)
/* calc the tangent of x */
/*!
* \brief Calculate the tangent of each element of a Matrix
*
* This function calculates the tangent of each element in a Matrix
*
* \param A The matrix whose tangents are of interest.
*
* \see sinh()
* \see tanh()
* \see sin()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(tan, ::tan)
/* calc the hyperbolic tangent of x */
/*!
* \brief Calculate the hyperbolic tangent of each element of a Matrix
*
* This function calculates the hyperbolic tangent of each element in a Matrix
*
* \param A The matrix whose hyperbolic tangents are of interest.
*
* \see sinh()
* \see tan()
* \see sin()
* \see cos()
* \see cosh()
* \see acos()
* \see acosh()
* \see asin()
* \see asinh()
* \see atan()
* \see atanh()
* \see atan2()
*/
SCYTHE_MATH_OP(tanh, ::tanh)
/* bessel function of the second kind of order 0*/
/*!
* \brief Compute the Bessel function of the second kind of order 0
*
* This function computes the Bessel function of the second kind of order 0
* for each element in the input matrix, A.
*
* \param A Matrix for which the Bessel function is of interest
*
* \see j0()
* \see j1()
* \see jn()
* \see y1()
* \see yn()
*/
SCYTHE_MATH_OP(y0, ::y0)
/* bessel function of the second kind of order 1*/
/*!
* \brief Compute the Bessel function of the second kind of order 1
*
* This function computes the Bessel function of the second kind of order 1
* for each element in the input matrix, A.
*
* \param A Matrix for which the Bessel function is of interest
*
* \see j0()
* \see j1()
* \see jn()
* \see y0()
* \see yn()
*/
SCYTHE_MATH_OP(y1, ::y1)
/* bessel function of the second kind of order n
* TODO: This definition causes the compiler to issue some warnings.
* Fix
*/
/*!
* \brief Compute the Bessel function of the second kind of order n
*
* This function computes the Bessel function of the second kind of order n
* for each element in the input matrix, A.
*
* \param n Order of the Bessel function
* \param A Matrix for which the Bessel function is of interest
*
* \see j0()
* \see j1()
* \see jn()
* \see y0()
* \see y1()
*/
SCYTHE_MATH_OP_2ARG(yn, std::ptr_fun(::yn))
} // end namespace scythe
#endif /* SCYTHE_MATH_H */
|