TABLE OF CONTENTS
- ABINIT/m_gaussian_quadrature
- m_hamiltonian/cdgqf
- m_hamiltonian/cgqf
- m_hamiltonian/ch_cap
- m_hamiltonian/ch_eqi
- m_hamiltonian/ch_to_digit
- m_hamiltonian/class_matrix
- m_hamiltonian/gaussian_quadrature_gegenbauer
- m_hamiltonian/gaussian_quadrature_legendre
- m_hamiltonian/get_frequencies_and_weights_legendre
- m_hamiltonian/imtqlx
- m_hamiltonian/parchk
- m_hamiltonian/r8_gamma_gq
- m_hamiltonian/r8mat_write
- m_hamiltonian/rule_write
- m_hamiltonian/s_to_i4
- m_hamiltonian/s_to_r8
- m_hamiltonian/scqf
- m_hamiltonian/sgqf
ABINIT/m_gaussian_quadrature [ Modules ]
NAME
m_gaussian_quadrature
FUNCTION
.
COPYRIGHT
Copyright (C) 2009-2022 ABINIT group (JLJ, BR, MC) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
SOURCE
16 #if defined HAVE_CONFIG_H 17 #include "config.h" 18 #endif 19 20 #include "abi_common.h" 21 22 module m_gaussian_quadrature 23 !---------------------------------------------------------------------------------------------------- 24 ! 25 ! This module contains routines to generate the weights and positions at which to evaluate 26 ! an integrand in order to produce a Gaussian quadrature of a given type. 27 ! 28 ! The code here was obtained from the internet (see header to subroutine cdgqf for further details). 29 ! The code was slightly modified to use doubles instead of floats in order to reach greater 30 ! accuracy. 31 ! 32 !---------------------------------------------------------------------------------------------------- 33 34 use m_abicore 35 use m_errors 36 37 use defs_basis, only : dp, tol20, std_out 38 use m_io_tools, only : open_file 39 40 implicit none 41 42 private
m_hamiltonian/cdgqf [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
cdgqf
FUNCTION
.
INPUTS
OUTPUT
SOURCE
332 subroutine cdgqf ( nt, gaussian_kind, alpha, beta, t, wts ) 333 334 !*****************************************************************************80 335 ! 336 !! CDGQF computes a Gauss quadrature formula with default A, B and simple knots. 337 ! 338 ! Discussion: 339 ! 340 ! This routine computes all the knots and weights of a Gauss quadrature 341 ! formula with a classical weight function with default values for A and B, 342 ! and only simple knots. 343 ! 344 ! There are no moments checks and no printing is done. 345 ! 346 ! Use routine EIQFS to evaluate a quadrature computed by CGQFS. 347 ! 348 ! Licensing: 349 ! 350 ! This code is distributed under the GNU LGPL license. 351 ! 352 ! Modified: 353 ! 354 ! 04 January 2010 355 ! 356 ! Author: 357 ! 358 ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. 359 ! FORTRAN90 version by John Burkardt. 360 ! 361 ! Reference: 362 ! 363 ! Sylvan Elhay, Jaroslav Kautsky, 364 ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of 365 ! Interpolatory Quadrature, 366 ! ACM Transactions on Mathematical Software, 367 ! Volume 13, Number 4, December 1987, pages 399-415. 368 ! 369 ! Parameters: 370 ! 371 ! Input, integer ( kind = 4 ) NT, the number of knots. 372 ! 373 ! Input, integer ( kind = 4 ) KIND, the rule. 374 ! 1, Legendre, (a,b) 1.0 375 ! 2, Chebyshev, (a,b) ((b-x)*(x-a))^(-0.5) 376 ! 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha 377 ! 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta 378 ! 5, Generalized Laguerre, (a,inf) (x-a)^alpha*exp(-b*(x-a)) 379 ! 6, Generalized Hermite, (-inf,inf) |x-a|^alpha*exp(-b*(x-a)^2) 380 ! 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha 381 ! 8, Rational, (a,inf) (x-a)^alpha*(x+b)^beta 382 ! 383 ! Input, real ( dp ) ALPHA, the value of Alpha, if needed. 384 ! 385 ! Input, real ( dp ) BETA, the value of Beta, if needed. 386 ! 387 ! Output, real ( dp ) T(NT), the knots. 388 ! 389 ! Output, real ( dp ) WTS(NT), the weights. 390 ! 391 integer ( kind = 4 ) nt 392 393 real ( dp ) aj(nt) 394 real ( dp ) alpha 395 real ( dp ) beta 396 real ( dp ) bj(nt) 397 integer ( kind = 4 ) gaussian_kind 398 real ( dp ) t(nt) 399 real ( dp ) wts(nt) 400 real ( dp ) zemu 401 402 ! ************************************************************************* 403 404 call parchk ( gaussian_kind, 2 * nt, alpha, beta ) 405 ! 406 ! Get the Jacobi matrix and zero-th moment. 407 ! 408 call class_matrix ( gaussian_kind, nt, alpha, beta, aj, bj, zemu ) 409 ! 410 ! Compute the knots and weights. 411 ! 412 call sgqf ( nt, aj, bj, zemu, t, wts ) 413 414 return 415 end subroutine cdgqf
m_hamiltonian/cgqf [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
cgqf
FUNCTION
.
INPUTS
OUTPUT
SOURCE
431 subroutine cgqf ( nt, gaussian_kind, alpha, beta, a, b, t, wts ) 432 433 !*****************************************************************************80 434 ! 435 !! CGQF computes knots and weights of a Gauss quadrature formula. 436 ! 437 ! Discussion: 438 ! 439 ! The user may specify the interval (A,B). 440 ! 441 ! Only simple knots are produced. 442 ! 443 ! Use routine EIQFS to evaluate this quadrature formula. 444 ! 445 ! Licensing: 446 ! 447 ! This code is distributed under the GNU LGPL license. 448 ! 449 ! Modified: 450 ! 451 ! 16 February 2010 452 ! 453 ! Author: 454 ! 455 ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. 456 ! FORTRAN90 version by John Burkardt. 457 ! 458 ! Reference: 459 ! 460 ! Sylvan Elhay, Jaroslav Kautsky, 461 ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of 462 ! Interpolatory Quadrature, 463 ! ACM Transactions on Mathematical Software, 464 ! Volume 13, Number 4, December 1987, pages 399-415. 465 ! 466 ! Parameters: 467 ! 468 ! Input, integer ( kind = 4 ) NT, the number of knots. 469 ! 470 ! Input, integer ( kind = 4 ) KIND, the rule. 471 ! 1, Legendre, (a,b) 1.0 472 ! 2, Chebyshev Type 1, (a,b) ((b-x)*(x-a))^-0.5) 473 ! 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha 474 ! 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta 475 ! 5, Generalized Laguerre, (a,+oo) (x-a)^alpha*exp(-b*(x-a)) 476 ! 6, Generalized Hermite, (-oo,+oo) |x-a|^alpha*exp(-b*(x-a)^2) 477 ! 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha 478 ! 8, Rational, (a,+oo) (x-a)^alpha*(x+b)^beta 479 ! 9, Chebyshev Type 2, (a,b) ((b-x)*(x-a))^(+0.5) 480 ! 481 ! Input, real ( dp ) ALPHA, the value of Alpha, if needed. 482 ! 483 ! Input, real ( dp ) BETA, the value of Beta, if needed. 484 ! 485 ! Input, real ( dp ) A, B, the interval endpoints, or 486 ! other parameters. 487 ! 488 ! Output, real ( dp ) T(NT), the knots. 489 ! 490 ! Output, real ( dp ) WTS(NT), the weights. 491 ! 492 integer ( kind = 4 ) nt 493 494 real ( dp ) a 495 real ( dp ) alpha 496 real ( dp ) b 497 real ( dp ) beta 498 integer ( kind = 4 ) i 499 integer ( kind = 4 ) gaussian_kind 500 integer ( kind = 4 ), allocatable :: mlt(:) 501 integer ( kind = 4 ), allocatable :: ndx(:) 502 real ( dp ) t(nt) 503 real ( dp ) wts(nt) 504 505 ! ************************************************************************* 506 507 ! 508 ! Compute the Gauss quadrature formula for default values of A and B. 509 ! 510 call cdgqf ( nt, gaussian_kind, alpha, beta, t, wts ) 511 ! 512 ! Prepare to scale the quadrature formula to other weight function with 513 ! valid A and B. 514 ! 515 ABI_MALLOC( mlt, (1:nt) ) 516 517 mlt(1:nt) = 1 518 519 ABI_MALLOC( ndx, (1:nt) ) 520 521 do i = 1, nt 522 ndx(i) = i 523 end do 524 525 call scqf ( nt, t, mlt, wts, nt, ndx, wts, t, gaussian_kind, alpha, beta, a, b ) 526 527 ABI_FREE( mlt ) 528 ABI_FREE( ndx ) 529 530 return 531 end subroutine cgqf
m_hamiltonian/ch_cap [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
ch_cap
FUNCTION
.
INPUTS
OUTPUT
SOURCE
547 subroutine ch_cap ( c ) 548 549 !*****************************************************************************80 550 ! 551 !! CH_CAP capitalizes a single character. 552 ! 553 ! Licensing: 554 ! 555 ! This code is distributed under the GNU LGPL license. 556 ! 557 ! Modified: 558 ! 559 ! 19 July 1998 560 ! 561 ! Author: 562 ! 563 ! John Burkardt 564 ! 565 ! Parameters: 566 ! 567 ! Input/output, character C, the character to capitalize. 568 ! 569 character c 570 integer ( kind = 4 ) itemp 571 572 ! ************************************************************************* 573 574 itemp = ichar ( c ) 575 576 if ( 97 <= itemp .and. itemp <= 122 ) then 577 c = char ( itemp - 32 ) 578 end if 579 580 return 581 end subroutine ch_cap
m_hamiltonian/ch_eqi [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
ch_eqi
FUNCTION
.
INPUTS
OUTPUT
SOURCE
599 function ch_eqi ( c1, c2 ) 600 601 !*****************************************************************************80 602 ! 603 !! CH_EQI is a case insensitive comparison of two characters for equality. 604 ! 605 ! Example: 606 ! 607 ! CH_EQI ( 'A', 'a' ) is .TRUE. 608 ! 609 ! Licensing: 610 ! 611 ! This code is distributed under the GNU LGPL license. 612 ! 613 ! Modified: 614 ! 615 ! 28 July 2000 616 ! 617 ! Author: 618 ! 619 ! John Burkardt 620 ! 621 ! Parameters: 622 ! 623 ! Input, character C1, C2, the characters to compare. 624 ! 625 ! Output, logical CH_EQI, the result of the comparison. 626 ! 627 logical ch_eqi 628 character c1 629 character c1_cap 630 character c2 631 character c2_cap 632 633 ! ************************************************************************* 634 635 c1_cap = c1 636 c2_cap = c2 637 638 call ch_cap ( c1_cap ) 639 call ch_cap ( c2_cap ) 640 641 if ( c1_cap == c2_cap ) then 642 ch_eqi = .true. 643 else 644 ch_eqi = .false. 645 end if 646 647 return 648 end function ch_eqi
m_hamiltonian/ch_to_digit [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
ch_to_digit
FUNCTION
.
INPUTS
OUTPUT
SOURCE
664 subroutine ch_to_digit ( c, digit ) 665 666 !*****************************************************************************80 667 ! 668 !! CH_TO_DIGIT returns the integer value of a base 10 digit. 669 ! 670 ! Example: 671 ! 672 ! C DIGIT 673 ! --- ----- 674 ! '0' 0 675 ! '1' 1 676 ! ... ... 677 ! '9' 9 678 ! ' ' 0 679 ! 'X' -1 680 ! 681 ! Licensing: 682 ! 683 ! This code is distributed under the GNU LGPL license. 684 ! 685 ! Modified: 686 ! 687 ! 04 August 1999 688 ! 689 ! Author: 690 ! 691 ! John Burkardt 692 ! 693 ! Parameters: 694 ! 695 ! Input, character C, the decimal digit, '0' through '9' or blank 696 ! are legal. 697 ! 698 ! Output, integer ( kind = 4 ) DIGIT, the corresponding integer value. 699 ! If C was 'illegal', then DIGIT is -1. 700 ! 701 character c 702 integer ( kind = 4 ) digit 703 704 ! ************************************************************************* 705 706 if ( lge ( c, '0' ) .and. lle ( c, '9' ) ) then 707 708 digit = ichar ( c ) - 48 709 710 else if ( c == ' ' ) then 711 712 digit = 0 713 714 else 715 716 digit = -1 717 718 end if 719 720 return 721 end subroutine ch_to_digit
m_hamiltonian/class_matrix [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
class_matrix
FUNCTION
.
INPUTS
OUTPUT
SOURCE
737 subroutine class_matrix ( gaussian_kind, m, alpha, beta, aj, bj, zemu ) 738 739 !*****************************************************************************80 740 ! 741 !! CLASS_MATRIX computes the Jacobi matrix for a quadrature rule. 742 ! 743 ! Discussion: 744 ! 745 ! This routine computes the diagonal AJ and sub-diagonal BJ 746 ! elements of the order M tridiagonal symmetric Jacobi matrix 747 ! associated with the polynomials orthogonal with respect to 748 ! the weight function specified by KIND. 749 ! 750 ! For weight functions 1-7, M elements are defined in BJ even 751 ! though only M-1 are needed. For weight function 8, BJ(M) is 752 ! set to zero. 753 ! 754 ! The zero-th moment of the weight function is returned in ZEMU. 755 ! 756 ! Licensing: 757 ! 758 ! This code is distributed under the GNU LGPL license. 759 ! 760 ! Modified: 761 ! 762 ! 27 December 2009 763 ! 764 ! Author: 765 ! 766 ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. 767 ! FORTRAN90 version by John Burkardt. 768 ! 769 ! Reference: 770 ! 771 ! Sylvan Elhay, Jaroslav Kautsky, 772 ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of 773 ! Interpolatory Quadrature, 774 ! ACM Transactions on Mathematical Software, 775 ! Volume 13, Number 4, December 1987, pages 399-415. 776 ! 777 ! Parameters: 778 ! 779 ! Input, integer ( kind = 4 ) KIND, the rule. 780 ! 1, Legendre, (a,b) 1.0 781 ! 2, Chebyshev, (a,b) ((b-x)*(x-a))^(-0.5) 782 ! 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha 783 ! 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta 784 ! 5, Generalized Laguerre, (a,inf) (x-a)^alpha*exp(-b*(x-a)) 785 ! 6, Generalized Hermite, (-inf,inf) |x-a|^alpha*exp(-b*(x-a)^2) 786 ! 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha 787 ! 8, Rational, (a,inf) (x-a)^alpha*(x+b)^beta 788 ! 789 ! Input, integer ( kind = 4 ) M, the order of the Jacobi matrix. 790 ! 791 ! Input, real ( dp ) ALPHA, the value of Alpha, if needed. 792 ! 793 ! Input, real ( dp ) BETA, the value of Beta, if needed. 794 ! 795 ! Output, real ( dp ) AJ(M), BJ(M), the diagonal and subdiagonal 796 ! of the Jacobi matrix. 797 ! 798 ! Output, real ( dp ) ZEMU, the zero-th moment. 799 ! 800 integer ( kind = 4 ) m 801 802 real ( dp ) a2b2 803 real ( dp ) ab 804 real ( dp ) aba 805 real ( dp ) abi 806 real ( dp ) abj 807 real ( dp ) abti 808 real ( dp ) aj(m) 809 real ( dp ) alpha 810 real ( dp ) apone 811 real ( dp ) beta 812 real ( dp ) bj(m) 813 integer ( kind = 4 ) i 814 integer ( kind = 4 ) gaussian_kind 815 real ( dp ), parameter :: pi = 3.14159265358979323846264338327950D+00 816 real ( dp ) temp 817 real ( dp ) temp2 818 real ( dp ) zemu 819 820 ! ************************************************************************* 821 822 temp = epsilon ( temp ) 823 824 call parchk ( gaussian_kind, 2 * m - 1, alpha, beta ) 825 826 temp2 = 0.5D+00 827 828 if ( 500.0D+00 * temp < abs ( ( r8_gamma_gq ( temp2 ) )**2 - pi ) ) then 829 ABI_ERROR('Gamma function does not match machine parameters.') 830 end if 831 832 if ( gaussian_kind == 1 ) then 833 834 ab = 0.0D+00 835 836 zemu = 2.0D+00 / ( ab + 1.0D+00 ) 837 838 aj(1:m) = 0.0D+00 839 840 do i = 1, m 841 abi = i + ab * mod ( i, 2 ) 842 abj = 2 * i + ab 843 bj(i) = abi * abi / ( abj * abj - 1.0D+00 ) 844 end do 845 bj(1:m) = sqrt ( bj(1:m) ) 846 847 else if ( gaussian_kind == 2 ) then 848 849 zemu = pi 850 851 aj(1:m) = 0.0D+00 852 853 bj(1) = sqrt ( 0.5D+00 ) 854 bj(2:m) = 0.5D+00 855 856 else if ( gaussian_kind == 3 ) then 857 858 ab = alpha * 2.0D+00 859 zemu = 2.0D+00**( ab + 1.0D+00 ) * r8_gamma_gq ( alpha + 1.0D+00 )**2 & 860 / r8_gamma_gq ( ab + 2.0D+00 ) 861 862 aj(1:m) = 0.0D+00 863 bj(1) = 1.0D+00 / ( 2.0D+00 * alpha + 3.0D+00 ) 864 do i = 2, m 865 bj(i) = i * ( i + ab ) / ( 4.0D+00 * ( i + alpha )**2 - 1.0D+00 ) 866 end do 867 bj(1:m) = sqrt ( bj(1:m) ) 868 869 else if ( gaussian_kind == 4 ) then 870 871 ab = alpha + beta 872 abi = 2.0D+00 + ab 873 zemu = 2.0D+00**( ab + 1.0D+00 ) * r8_gamma_gq ( alpha + 1.0D+00 ) & 874 * r8_gamma_gq ( beta + 1.0D+00 ) / r8_gamma_gq ( abi ) 875 aj(1) = ( beta - alpha ) / abi 876 bj(1) = 4.0D+00 * ( 1.0 + alpha ) * ( 1.0D+00 + beta ) & 877 / ( ( abi + 1.0D+00 ) * abi * abi ) 878 a2b2 = beta * beta - alpha * alpha 879 880 do i = 2, m 881 abi = 2.0D+00 * i + ab 882 aj(i) = a2b2 / ( ( abi - 2.0D+00 ) * abi ) 883 abi = abi**2 884 bj(i) = 4.0D+00 * i * ( i + alpha ) * ( i + beta ) * ( i + ab ) & 885 / ( ( abi - 1.0D+00 ) * abi ) 886 end do 887 bj(1:m) = sqrt ( bj(1:m) ) 888 889 else if ( gaussian_kind == 5 ) then 890 891 zemu = r8_gamma_gq ( alpha + 1.0D+00 ) 892 893 do i = 1, m 894 aj(i) = 2.0D+00 * i - 1.0D+00 + alpha 895 bj(i) = i * ( i + alpha ) 896 end do 897 bj(1:m) = sqrt ( bj(1:m) ) 898 899 else if ( gaussian_kind == 6 ) then 900 901 zemu = r8_gamma_gq ( ( alpha + 1.0D+00 ) / 2.0D+00 ) 902 903 aj(1:m) = 0.0D+00 904 905 do i = 1, m 906 bj(i) = ( i + alpha * mod ( i, 2 ) ) / 2.0D+00 907 end do 908 bj(1:m) = sqrt ( bj(1:m) ) 909 910 else if ( gaussian_kind == 7 ) then 911 912 ab = alpha 913 zemu = 2.0D+00 / ( ab + 1.0D+00 ) 914 915 aj(1:m) = 0.0D+00 916 917 do i = 1, m 918 abi = i + ab * mod ( i, 2 ) 919 abj = 2 * i + ab 920 bj(i) = abi * abi / ( abj * abj - 1.0D+00 ) 921 end do 922 bj(1:m) = sqrt ( bj(1:m) ) 923 924 else if ( gaussian_kind == 8 ) then 925 926 ab = alpha + beta 927 zemu = r8_gamma_gq ( alpha + 1.0D+00 ) * r8_gamma_gq ( - ( ab + 1.0D+00 ) ) & 928 / r8_gamma_gq ( - beta ) 929 apone = alpha + 1.0D+00 930 aba = ab * apone 931 aj(1) = - apone / ( ab + 2.0D+00 ) 932 bj(1) = - aj(1) * ( beta + 1.0D+00 ) / ( ab + 2.0D+00 ) / ( ab + 3.0D+00 ) 933 do i = 2, m 934 abti = ab + 2.0D+00 * i 935 aj(i) = aba + 2.0D+00 * ( ab + i ) * ( i - 1 ) 936 aj(i) = - aj(i) / abti / ( abti - 2.0D+00 ) 937 end do 938 939 do i = 2, m - 1 940 abti = ab + 2.0D+00 * i 941 bj(i) = i * ( alpha + i ) / ( abti - 1.0D+00 ) * ( beta + i ) & 942 / ( abti**2 ) * ( ab + i ) / ( abti + 1.0D+00 ) 943 end do 944 945 bj(m) = 0.0D+00 946 bj(1:m) = sqrt ( bj(1:m) ) 947 948 end if 949 950 return 951 end subroutine class_matrix
m_hamiltonian/gaussian_quadrature_gegenbauer [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
gaussian_quadrature_gegenbauer
FUNCTION
.
INPUTS
OUTPUT
SOURCE
67 subroutine gaussian_quadrature_gegenbauer(integrand,integral,order) 68 !---------------------------------------------------------------------------------------------------- 69 ! 70 ! This subroutine performs a gaussian quadrature of given order of the input function integrand 71 ! and stores the result in the variable integral. 72 ! 73 ! The integrand function is assumed to be defined on the range from 0 to infinity, and it is 74 ! assumed to to behave as x^2 near x=0, and 1/x^4 as x-> infinity. 75 ! 76 ! The Gegenbauer integration method is of the form: 77 ! 78 ! int_{-1}^{1} dt (1-t)^2(1+t)^2 f(t) = sum_{i} w_i f(t_i) 79 ! 80 ! 81 ! In order to perform such an integral, make a change of variable x = (1-t)/(1+t); this implies 82 ! 83 ! int_0^infty dx I(x) = int_{-1}^{1} dt g(t), g(t) = 2/(1+t)^2 I((1-t)/(1+t)) 84 ! 85 ! It is easy to show that g(1-delta) ~ delta^2 and g(-1+delta) ~ delta^2 for delta small. 86 ! Thus, we can define 87 ! 88 ! int_{-1}^{1} dt g(t) = int_{-1}^{1} dt (1+t)^2(1-t)^2 f(t), f(t) = g(t)/(1-t^2)^2 89 ! 90 ! 91 !---------------------------------------------------------------------------------------------------- 92 93 94 ! input and output variables 95 integer, intent(in) :: order 96 97 interface 98 function integrand(x) 99 integer, parameter :: dp=kind(1.0d0) 100 real(dp) integrand 101 real(dp),intent(in) :: x 102 end function integrand 103 end interface 104 105 real(dp), intent(out) :: integral 106 107 108 real(dp) :: t(order), w(order) 109 110 111 integer :: i 112 113 integer :: gaussian_kind 114 real(dp) :: alpha, beta, a, b 115 real(dp) :: x, y, f 116 117 ! ************************************************************************* 118 119 120 ! Set the arguments and weights 121 gaussian_kind = 3 122 alpha= 2.0_dp 123 beta = 2.0_dp 124 a = -1.0_dp 125 b = 1.0_dp 126 call cgqf ( order, gaussian_kind, alpha, beta, a, b, t, w ) 127 128 129 130 ! call the integrant the right number of times 131 132 integral = 0.0_dp 133 134 do i = 1, order 135 x = (1.0_dp-t(i))/(1.0_dp+t(i)) 136 137 y = integrand(x) 138 139 f = 2.0_dp/(1.0_dp+t(i))**2/(1.0_dp-t(i)**2)**2*y 140 141 integral = integral + w(i)*f 142 143 end do ! i 144 145 end subroutine gaussian_quadrature_gegenbauer
m_hamiltonian/gaussian_quadrature_legendre [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
gaussian_quadrature_legendre
FUNCTION
.
INPUTS
OUTPUT
SOURCE
161 subroutine gaussian_quadrature_legendre(integrand,integral,order) 162 !---------------------------------------------------------------------------------------------------- 163 ! 164 ! This subroutine performs a gaussian quadrature of given order of the input function integrand 165 ! and stores the result in the variable integral. 166 ! 167 ! The integrand function is assumed to be defined on the range from 0 to infinity, and it is 168 ! assumed to to behave as x^2 near x=0, and 1/x^4 as x-> infinity. 169 ! 170 ! The Legendre integration method is of the form: 171 ! 172 ! int_{-1}^{1} dt f(t) = sum_{i} w_i f(t_i) 173 ! 174 ! 175 ! In order to perform such an integral, make a change of variable x = (1-t)/(1+t); this implies 176 ! 177 ! int_0^infty dx I(x) = int_{-1}^{1} dt f(t), f(t) = 2/(1+t)^2 I((1-t)/(1+t)) 178 ! 179 ! 180 ! 181 !---------------------------------------------------------------------------------------------------- 182 183 184 ! input and output variables 185 integer, intent(in) :: order 186 187 interface 188 function integrand(x) 189 integer, parameter :: dp=kind(1.0d0) 190 real(dp) integrand 191 real(dp),intent(in) :: x 192 end function integrand 193 end interface 194 195 real(dp), intent(out) :: integral 196 197 198 real(dp) :: t(order), w(order) 199 200 201 integer :: i 202 203 integer :: gaussian_kind 204 real(dp) :: alpha, beta, a, b 205 real(dp) :: x, y, f 206 207 ! ************************************************************************* 208 209 210 ! Set the arguments and weights 211 gaussian_kind = 1 212 alpha= 0.0_dp 213 beta = 0.0_dp 214 a = -1.0_dp 215 b = 1.0_dp 216 call cgqf ( order, gaussian_kind, alpha, beta, a, b, t, w ) 217 218 219 220 ! call the integrant the right number of times 221 222 integral = 0.0_dp 223 224 do i = 1, order 225 x = (1.0_dp-t(i))/(1.0_dp+t(i)) 226 227 y = integrand(x) 228 229 f = 2.0_dp/(1.0_dp+t(i))**2*y 230 231 integral = integral + w(i)*f 232 233 end do ! i 234 235 end subroutine gaussian_quadrature_legendre
m_hamiltonian/get_frequencies_and_weights_legendre [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
get_frequencies_and_weights_legendre
FUNCTION
.
INPUTS
OUTPUT
SOURCE
252 subroutine get_frequencies_and_weights_legendre(order,frequencies,weights) 253 !---------------------------------------------------------------------------------------------------- 254 ! 255 ! This subroutine generates the frequencies and weights to integrate using legendre Gaussian 256 ! quadrature for a given order. 257 ! 258 ! The integrand function is assumed to be defined on the range from 0 to infinity, and it is 259 ! assumed to to behave as x^2 near x=0, and 1/x^4 as x-> infinity. 260 ! 261 ! The Legendre integration method is of the form: 262 ! 263 ! int_{-1}^{1} dt f(t) = sum_{i} w_i f(t_i) 264 ! 265 ! In order to perform such an integral, make a change of variable 266 ! 267 ! omega = (1-t)/(1+t) 268 ! 269 ! and thus 270 ! int_0^infty domega I(omega) = int_{-1}^{1} dt f(t), f(t) = 2/(1+t)^2 I((1-t)/(1+t)) 271 ! 272 ! and so the weights are defined as weights(i) = 2/(1+ti)^2*wi. 273 ! 274 ! 275 !---------------------------------------------------------------------------------------------------- 276 ! input and output variables 277 integer, intent(in) :: order 278 real(dp), intent(out) :: frequencies(order) 279 real(dp), intent(out) :: weights(order) 280 281 282 real(dp) :: t(order), w(order) 283 284 285 integer :: i 286 287 integer :: gaussian_kind 288 real(dp) :: alpha, beta, a, b 289 real(dp) :: x, f 290 291 ! ************************************************************************* 292 293 294 ! Set the arguments and weights 295 gaussian_kind = 1 296 alpha = 0.0_dp 297 beta = 0.0_dp 298 a = -1.0_dp 299 b = 1.0_dp 300 301 call cgqf ( order, gaussian_kind, alpha, beta, a, b, t, w ) 302 303 304 ! call the integrant the right number of times 305 306 do i = 1, order 307 308 x = (1.0_dp-t(i))/(1.0_dp+t(i)) 309 f = 2.0_dp/(1.0_dp+t(i))**2 310 311 frequencies(i) = x 312 weights(i) = w(i)*f 313 314 end do ! i 315 316 end subroutine get_frequencies_and_weights_legendre
m_hamiltonian/imtqlx [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
imtqlx
FUNCTION
.
INPUTS
OUTPUT
SOURCE
967 subroutine imtqlx ( n, d, e, z ) 968 969 !*****************************************************************************80 970 ! 971 !! IMTQLX diagonalizes a symmetric tridiagonal matrix. 972 ! 973 ! Discussion: 974 ! 975 ! This routine is a slightly modified version of the EISPACK routine to 976 ! perform the implicit QL algorithm on a symmetric tridiagonal matrix. 977 ! 978 ! The authors thank the authors of EISPACK for permission to use this 979 ! routine. 980 ! 981 ! It has been modified to produce the product Q' * Z, where Z is an input 982 ! vector and Q is the orthogonal matrix diagonalizing the input matrix. 983 ! The changes consist (essentially) of applying the orthogonal 984 ! transformations directly to Z as they are generated. 985 ! 986 ! Licensing: 987 ! 988 ! This code is distributed under the GNU LGPL license. 989 ! 990 ! Modified: 991 ! 992 ! 27 December 2009 993 ! 994 ! Author: 995 ! 996 ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. 997 ! FORTRAN90 version by John Burkardt. 998 ! 999 ! Reference: 1000 ! 1001 ! Sylvan Elhay, Jaroslav Kautsky, 1002 ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of 1003 ! Interpolatory Quadrature, 1004 ! ACM Transactions on Mathematical Software, 1005 ! Volume 13, Number 4, December 1987, pages 399-415. 1006 ! 1007 ! Roger Martin, James Wilkinson, 1008 ! The Implicit QL Algorithm, 1009 ! Numerische Mathematik, 1010 ! Volume 12, Number 5, December 1968, pages 377-383. 1011 ! 1012 ! Parameters: 1013 ! 1014 ! Input, integer ( kind = 4 ) N, the order of the matrix. 1015 ! 1016 ! Input/output, real ( dp ) D(N), the diagonal entries of the matrix. 1017 ! On output, the information in D has been overwritten. 1018 ! 1019 ! Input/output, real ( dp ) E(N), the subdiagonal entries of the 1020 ! matrix, in entries E(1) through E(N-1). On output, the information in 1021 ! E has been overwritten. 1022 ! 1023 ! Input/output, real ( dp ) Z(N). On input, a vector. On output, 1024 ! the value of Q' * Z, where Q is the matrix that diagonalizes the 1025 ! input symmetric tridiagonal matrix. 1026 ! 1027 integer ( kind = 4 ) n 1028 1029 real ( dp ) b 1030 real ( dp ) c 1031 real ( dp ) d(n) 1032 real ( dp ) e(n) 1033 real ( dp ) f 1034 real ( dp ) g 1035 integer ( kind = 4 ) i 1036 integer ( kind = 4 ) ii 1037 integer ( kind = 4 ), parameter :: itn = 30 1038 integer ( kind = 4 ) j 1039 integer ( kind = 4 ) k 1040 integer ( kind = 4 ) l 1041 integer ( kind = 4 ) m 1042 integer ( kind = 4 ) mml 1043 real ( dp ) p 1044 real ( dp ) prec 1045 real ( dp ) r 1046 real ( dp ) s 1047 real ( dp ) z(n) 1048 1049 ! ************************************************************************* 1050 1051 prec = epsilon ( prec ) 1052 1053 if ( n == 1 ) then 1054 return 1055 end if 1056 1057 e(n) = 0.0D+00 1058 1059 do l = 1, n 1060 1061 j = 0 1062 1063 do 1064 1065 do m = l, n 1066 1067 if ( m == n ) then 1068 exit 1069 end if 1070 1071 if ( abs ( e(m) ) <= prec * ( abs ( d(m) ) + abs ( d(m+1) ) ) ) then 1072 exit 1073 end if 1074 1075 end do 1076 1077 p = d(l) 1078 1079 if ( m == l ) then 1080 exit 1081 end if 1082 1083 if ( itn <= j ) then 1084 write ( std_out, '(a)' ) ' ' 1085 write ( std_out, '(a)' ) 'IMTQLX - Fatal error!' 1086 write ( std_out, '(a)' ) ' Iteration limit exceeded.' 1087 write ( std_out, '(a,i8)' ) ' J = ', j 1088 write ( std_out, '(a,i8)' ) ' L = ', l 1089 write ( std_out, '(a,i8)' ) ' M = ', m 1090 write ( std_out, '(a,i8)' ) ' N = ', n 1091 ABI_ERROR("Aborting now") 1092 end if 1093 1094 j = j + 1 1095 g = ( d(l+1) - p ) / ( 2.0D+00 * e(l) ) 1096 r = sqrt ( g * g + 1.0D+00 ) 1097 g = d(m) - p + e(l) / ( g + sign ( r, g ) ) 1098 s = 1.0D+00 1099 c = 1.0D+00 1100 p = 0.0D+00 1101 mml = m - l 1102 1103 do ii = 1, mml 1104 1105 i = m - ii 1106 f = s * e(i) 1107 b = c * e(i) 1108 1109 if ( abs ( g ) <= abs ( f ) ) then 1110 c = g / f 1111 r = sqrt ( c * c + 1.0D+00 ) 1112 e(i+1) = f * r 1113 s = 1.0D+00 / r 1114 c = c * s 1115 else 1116 s = f / g 1117 r = sqrt ( s * s + 1.0D+00 ) 1118 e(i+1) = g * r 1119 c = 1.0D+00 / r 1120 s = s * c 1121 end if 1122 1123 g = d(i+1) - p 1124 r = ( d(i) - g ) * s + 2.0D+00 * c * b 1125 p = s * r 1126 d(i+1) = g + p 1127 g = c * r - b 1128 f = z(i+1) 1129 z(i+1) = s * z(i) + c * f 1130 z(i) = c * z(i) - s * f 1131 1132 end do 1133 1134 d(l) = d(l) - p 1135 e(l) = g 1136 e(m) = 0.0D+00 1137 1138 end do 1139 1140 end do 1141 ! 1142 ! Sorting. 1143 ! 1144 do ii = 2, n 1145 1146 i = ii - 1 1147 k = i 1148 p = d(i) 1149 1150 do j = ii, n 1151 if ( d(j) < p ) then 1152 k = j 1153 p = d(j) 1154 end if 1155 end do 1156 1157 if ( k /= i ) then 1158 d(k) = d(i) 1159 d(i) = p 1160 p = z(i) 1161 z(i) = z(k) 1162 z(k) = p 1163 end if 1164 1165 end do 1166 1167 return 1168 end subroutine imtqlx
m_hamiltonian/parchk [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
parchk
FUNCTION
.
INPUTS
OUTPUT
SOURCE
1184 subroutine parchk ( gaussian_kind, m, alpha, beta ) 1185 1186 !*****************************************************************************80 1187 ! 1188 !! PARCHK checks parameters ALPHA and BETA for classical weight functions. 1189 ! 1190 ! Licensing: 1191 ! 1192 ! This code is distributed under the GNU LGPL license. 1193 ! 1194 ! Modified: 1195 ! 1196 ! 27 December 2009 1197 ! 1198 ! Author: 1199 ! 1200 ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. 1201 ! FORTRAN90 version by John Burkardt. 1202 ! 1203 ! Reference: 1204 ! 1205 ! Sylvan Elhay, Jaroslav Kautsky, 1206 ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of 1207 ! Interpolatory Quadrature, 1208 ! ACM Transactions on Mathematical Software, 1209 ! Volume 13, Number 4, December 1987, pages 399-415. 1210 ! 1211 ! Parameters: 1212 ! 1213 ! Input, integer ( kind = 4 ) KIND, the rule. 1214 ! 1, Legendre, (a,b) 1.0 1215 ! 2, Chebyshev, (a,b) ((b-x)*(x-a))^(-0.5) 1216 ! 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha 1217 ! 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta 1218 ! 5, Generalized Laguerre, (a,inf) (x-a)^alpha*exp(-b*(x-a)) 1219 ! 6, Generalized Hermite, (-inf,inf) |x-a|^alpha*exp(-b*(x-a)^2) 1220 ! 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha 1221 ! 8, Rational, (a,inf) (x-a)^alpha*(x+b)^beta 1222 ! 1223 ! Input, integer ( kind = 4 ) M, the order of the highest moment to 1224 ! be calculated. This value is only needed when KIND = 8. 1225 ! 1226 ! Input, real ( dp ) ALPHA, BETA, the parameters, if required 1227 ! by the value of KIND. 1228 ! 1229 real ( dp ) alpha 1230 real ( dp ) beta 1231 integer ( kind = 4 ) gaussian_kind 1232 integer ( kind = 4 ) m 1233 real ( dp ) tmp 1234 1235 ! ************************************************************************* 1236 1237 if ( gaussian_kind <= 0 ) then 1238 ABI_ERROR('PARCHK - Fatal error: KIND <= 0.') 1239 end if 1240 ! 1241 ! Check ALPHA for Gegenbauer, Jacobi, Laguerre, Hermite, Exponential. 1242 ! 1243 if ( 3 <= gaussian_kind .and. alpha <= -1.0D+00 ) then 1244 ABI_ERROR('PARCHK - Fatal error! 3 <= KIND and ALPHA <= -1.') 1245 end if 1246 ! 1247 ! Check BETA for Jacobi. 1248 ! 1249 if ( gaussian_kind == 4 .and. beta <= -1.0D+00 ) then 1250 ABI_ERROR('PARCHK - Fatal error: KIND == 4 and BETA <= -1.0.') 1251 end if 1252 ! 1253 ! Check ALPHA and BETA for rational. 1254 ! 1255 if ( gaussian_kind == 8 ) then 1256 tmp = alpha + beta + m + 1.0D+00 1257 if ( 0.0D+00 <= tmp .or. tmp <= beta ) then 1258 ABI_ERROR('PARCHK - Fatal error! KIND == 8 but condition on ALPHA and BETA fails.') 1259 end if 1260 end if 1261 1262 return 1263 end subroutine parchk
m_hamiltonian/r8_gamma_gq [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
r8_gamma_gq
FUNCTION
.
INPUTS
OUTPUT
SOURCE
1282 function r8_gamma_gq ( x ) 1283 1284 !*****************************************************************************80 1285 ! 1286 !! R8_GAMMA evaluates Gamma(X) for a real argument. 1287 ! 1288 ! Discussion: 1289 ! 1290 ! This routine calculates the gamma function for a real argument X. 1291 ! 1292 ! Computation is based on an algorithm outlined in reference 1. 1293 ! The program uses rational functions that approximate the gamma 1294 ! function to at least 20 significant decimal digits. Coefficients 1295 ! for the approximation over the interval (1,2) are unpublished. 1296 ! Those for the approximation for 12 <= X are from reference 2. 1297 ! 1298 ! Licensing: 1299 ! 1300 ! This code is distributed under the GNU LGPL license. 1301 ! 1302 ! Modified: 1303 ! 1304 ! 11 February 2008 1305 ! 1306 ! Author: 1307 ! 1308 ! Original FORTRAN77 version by William Cody, Laura Stoltz. 1309 ! FORTRAN90 version by John Burkardt. 1310 ! 1311 ! Reference: 1312 ! 1313 ! William Cody, 1314 ! An Overview of Software Development for Special Functions, 1315 ! in Numerical Analysis Dundee, 1975, 1316 ! edited by GA Watson, 1317 ! Lecture Notes in Mathematics 506, 1318 ! Springer, 1976. 1319 ! 1320 ! John Hart, Ward Cheney, Charles Lawson, Hans Maehly, 1321 ! Charles Mesztenyi, John Rice, Henry Thatcher, 1322 ! Christoph Witzgall, 1323 ! Computer Approximations, 1324 ! Wiley, 1968, 1325 ! LC: QA297.C64. 1326 ! 1327 ! Parameters: 1328 ! 1329 ! Input, real ( dp ) X, the argument of the function. 1330 ! 1331 ! Output, real ( dp ) R8_GAMMA, the value of the function. 1332 ! 1333 real ( dp ), dimension ( 7 ) :: c = (/ & 1334 -1.910444077728D-03, & 1335 8.4171387781295D-04, & 1336 -5.952379913043012D-04, & 1337 7.93650793500350248D-04, & 1338 -2.777777777777681622553D-03, & 1339 8.333333333333333331554247D-02, & 1340 5.7083835261D-03 /) 1341 real ( dp ), parameter :: eps = 2.22D-16 1342 real ( dp ) fact 1343 integer ( kind = 4 ) i 1344 integer ( kind = 4 ) n 1345 real ( dp ), dimension ( 8 ) :: p = (/ & 1346 -1.71618513886549492533811D+00, & 1347 2.47656508055759199108314D+01, & 1348 -3.79804256470945635097577D+02, & 1349 6.29331155312818442661052D+02, & 1350 8.66966202790413211295064D+02, & 1351 -3.14512729688483675254357D+04, & 1352 -3.61444134186911729807069D+04, & 1353 6.64561438202405440627855D+04 /) 1354 logical parity 1355 real ( dp ), parameter :: pi = 3.1415926535897932384626434D+00 1356 real ( dp ), dimension ( 8 ) :: q = (/ & 1357 -3.08402300119738975254353D+01, & 1358 3.15350626979604161529144D+02, & 1359 -1.01515636749021914166146D+03, & 1360 -3.10777167157231109440444D+03, & 1361 2.25381184209801510330112D+04, & 1362 4.75584627752788110767815D+03, & 1363 -1.34659959864969306392456D+05, & 1364 -1.15132259675553483497211D+05 /) 1365 real ( dp ) r8_gamma_gq 1366 real ( dp ) res 1367 real ( dp ), parameter :: sqrtpi = 0.9189385332046727417803297D+00 1368 real ( dp ) sum 1369 real ( dp ) x 1370 real ( dp ), parameter :: xbig = 171.624D+00 1371 real ( dp ) xden 1372 real ( dp ), parameter :: xinf = 1.0D+30 1373 real ( dp ), parameter :: xminin = 2.23D-308 1374 real ( dp ) xnum 1375 real ( dp ) y 1376 real ( dp ) y1 1377 real ( dp ) ysq 1378 real ( dp ) z 1379 1380 ! ************************************************************************* 1381 1382 parity = .false. 1383 fact = 1.0D+00 1384 n = 0 1385 y = x 1386 ! 1387 ! Argument is negative. 1388 ! 1389 if ( y <= 0.0D+00 ) then 1390 1391 y = - x 1392 y1 = aint ( y ) 1393 res = y - y1 1394 1395 if (abs(res) > tol20) then !if ( res /= 0.0D+00 ) then 1396 1397 if ( abs(y1 - aint ( y1 * 0.5D+00 ) * 2.0D+00) > tol20 ) then ! if ( y1 /= aint ( y1 * 0.5D+00 ) * 2.0D+00 ) then 1398 parity = .true. 1399 end if 1400 1401 fact = - pi / sin ( pi * res ) 1402 y = y + 1.0D+00 1403 1404 else 1405 1406 res = xinf 1407 r8_gamma_gq = res 1408 return 1409 1410 end if 1411 1412 end if 1413 ! 1414 ! Argument is positive. 1415 ! 1416 if ( y < eps ) then 1417 ! 1418 ! Argument < EPS. 1419 ! 1420 if ( xminin <= y ) then 1421 res = 1.0D+00 / y 1422 else 1423 res = xinf 1424 r8_gamma_gq = res 1425 return 1426 end if 1427 1428 else if ( y < 12.0D+00 ) then 1429 1430 y1 = y 1431 ! 1432 ! 0.0 < argument < 1.0. 1433 ! 1434 if ( y < 1.0D+00 ) then 1435 1436 z = y 1437 y = y + 1.0D+00 1438 ! 1439 ! 1.0 < argument < 12.0. 1440 ! Reduce argument if necessary. 1441 ! 1442 else 1443 1444 n = int ( y ) - 1 1445 y = y - real ( n, dp ) 1446 z = y - 1.0D+00 1447 1448 end if 1449 ! 1450 ! Evaluate approximation for 1.0 < argument < 2.0. 1451 ! 1452 xnum = 0.0D+00 1453 xden = 1.0D+00 1454 do i = 1, 8 1455 xnum = ( xnum + p(i) ) * z 1456 xden = xden * z + q(i) 1457 end do 1458 1459 res = xnum / xden + 1.0D+00 1460 ! 1461 ! Adjust result for case 0.0 < argument < 1.0. 1462 ! 1463 if ( y1 < y ) then 1464 1465 res = res / y1 1466 ! 1467 ! Adjust result for case 2.0 < argument < 12.0. 1468 ! 1469 else if ( y < y1 ) then 1470 1471 do i = 1, n 1472 res = res * y 1473 y = y + 1.0D+00 1474 end do 1475 1476 end if 1477 1478 else 1479 ! 1480 ! Evaluate for 12.0 <= argument. 1481 ! 1482 if ( y <= xbig ) then 1483 1484 ysq = y * y 1485 sum = c(7) 1486 do i = 1, 6 1487 sum = sum / ysq + c(i) 1488 end do 1489 sum = sum / y - y + sqrtpi 1490 sum = sum + ( y - 0.5D+00 ) * log ( y ) 1491 res = exp ( sum ) 1492 1493 else 1494 1495 res = xinf 1496 r8_gamma_gq = res 1497 return 1498 1499 end if 1500 1501 end if 1502 ! 1503 ! Final adjustments and return. 1504 ! 1505 if ( parity ) then 1506 res = - res 1507 end if 1508 1509 if ( abs(fact - 1.0D+00) > tol20 ) then !if ( fact /= 1.0D+00 ) then 1510 res = fact / res 1511 end if 1512 1513 r8_gamma_gq = res 1514 1515 return 1516 end function r8_gamma_gq
m_hamiltonian/r8mat_write [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
r8mat_write
FUNCTION
.
INPUTS
OUTPUT
SOURCE
1532 subroutine r8mat_write ( output_filename, m, n, table ) 1533 1534 !*****************************************************************************80 1535 ! 1536 !! R8MAT_WRITE writes an R8MAT file. 1537 ! 1538 ! Licensing: 1539 ! 1540 ! This code is distributed under the GNU LGPL license. 1541 ! 1542 ! Modified: 1543 ! 1544 ! 31 May 2009 1545 ! 1546 ! Author: 1547 ! 1548 ! John Burkardt 1549 ! 1550 ! Parameters: 1551 ! 1552 ! Input, character ( len = * ) OUTPUT_FILENAME, the output file name. 1553 ! 1554 ! Input, integer ( kind = 4 ) M, the spatial dimension. 1555 ! 1556 ! Input, integer ( kind = 4 ) N, the number of points. 1557 ! 1558 ! Input, real ( dp ) TABLE(M,N), the table data. 1559 ! 1560 integer ( kind = 4 ) m 1561 integer ( kind = 4 ) n 1562 1563 integer ( kind = 4 ) j 1564 character ( len = * ) output_filename 1565 character(len=500) :: msg 1566 integer ( kind = 4 ) output_unit 1567 character ( len = 30 ) string 1568 real ( dp ) table(m,n) 1569 1570 ! ************************************************************************* 1571 1572 ! 1573 ! Open the file. 1574 ! 1575 if (open_file(output_filename, msg, newunit=output_unit, status='replace') /= 0) then 1576 ABI_ERROR(msg) 1577 end if 1578 1579 ! 1580 ! Create a format string. 1581 ! 1582 ! For less precision in the output file, try: 1583 ! 1584 ! '(', m, 'g', 14, '.', 6, ')' 1585 ! 1586 if ( 0 < m .and. 0 < n ) then 1587 1588 write ( string, '(a1,i8,a1,i8,a1,i8,a1)' ) '(', m, 'g', 24, '.', 16, ')' 1589 ! 1590 ! Write the data. 1591 ! 1592 do j = 1, n 1593 write ( output_unit, string ) table(1:m,j) 1594 end do 1595 1596 end if 1597 ! 1598 ! Close the file. 1599 ! 1600 close ( unit = output_unit ) 1601 1602 end subroutine r8mat_write
m_hamiltonian/rule_write [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
rule_write
FUNCTION
.
INPUTS
OUTPUT
SOURCE
1619 subroutine rule_write ( order, x, w, r, filename ) 1620 1621 !*****************************************************************************80 1622 ! 1623 !! RULE_WRITE writes a quadrature rule to a file. 1624 ! 1625 ! Licensing: 1626 ! 1627 ! This code is distributed under the GNU LGPL license. 1628 ! 1629 ! Modified: 1630 ! 1631 ! 18 February 2010 1632 ! 1633 ! Author: 1634 ! 1635 ! John Burkardt 1636 ! 1637 ! Parameters: 1638 ! 1639 ! Input, integer ( kind = 4 ) ORDER, the order of the rule. 1640 ! 1641 ! Input, real ( dp ) X(ORDER), the abscissas. 1642 ! 1643 ! Input, real ( dp ) W(ORDER), the weights. 1644 ! 1645 ! Input, real ( dp ) R(2), defines the region. 1646 ! 1647 ! Input, character ( len = * ) FILENAME, specifies the output. 1648 ! 'filename_w.txt', 'filename_x.txt', 'filename_r.txt' defining weights, 1649 ! abscissas, and region. 1650 ! 1651 integer ( kind = 4 ) order 1652 1653 character ( len = * ) filename 1654 character ( len = 255 ) filename_r 1655 character ( len = 255 ) filename_w 1656 character ( len = 255 ) filename_x 1657 real ( dp ) r(2) 1658 real ( dp ) w(order) 1659 real ( dp ) x(order) 1660 1661 ! ************************************************************************* 1662 1663 filename_w = trim ( filename ) // '_w.txt' 1664 filename_x = trim ( filename ) // '_x.txt' 1665 filename_r = trim ( filename ) // '_r.txt' 1666 1667 write ( std_out, '(a)' ) ' ' 1668 write ( std_out, '(a)' ) ' Creating quadrature files.' 1669 write ( std_out, '(a)' ) ' ' 1670 write ( std_out, '(a)' ) ' "Root" file name is "' // trim ( filename ) // '".' 1671 write ( std_out, '(a)' ) ' ' 1672 write ( std_out, '(a)' ) ' Weight file will be "' // trim ( filename_w ) // '".' 1673 write ( std_out, '(a)' ) ' Abscissa file will be "' // trim ( filename_x ) // '".' 1674 write ( std_out, '(a)' ) ' Region file will be "' // trim ( filename_r ) // '".' 1675 1676 call r8mat_write ( filename_w, 1, order, w ) 1677 call r8mat_write ( filename_x, 1, order, x ) 1678 call r8mat_write ( filename_r, 1, 2, r ) 1679 1680 return 1681 end subroutine rule_write
m_hamiltonian/s_to_i4 [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
s_to_i4
FUNCTION
.
INPUTS
OUTPUT
SOURCE
1698 subroutine s_to_i4 ( s, ival, ierror, length ) 1699 1700 !*****************************************************************************80 1701 ! 1702 !! S_TO_I4 reads an I4 from a string. 1703 ! 1704 ! Licensing: 1705 ! 1706 ! This code is distributed under the GNU LGPL license. 1707 ! 1708 ! Modified: 1709 ! 1710 ! 28 June 2000 1711 ! 1712 ! Author: 1713 ! 1714 ! John Burkardt 1715 ! 1716 ! Parameters: 1717 ! 1718 ! Input, character ( len = * ) S, a string to be examined. 1719 ! 1720 ! Output, integer ( kind = 4 ) IVAL, the integer value read from the string. 1721 ! If the string is blank, then IVAL will be returned 0. 1722 ! 1723 ! Output, integer ( kind = 4 ) IERROR, an error flag. 1724 ! 0, no error. 1725 ! 1, an error occurred. 1726 ! 1727 ! Output, integer ( kind = 4 ) LENGTH, the number of characters of S 1728 ! used to make IVAL. 1729 ! 1730 character c 1731 integer ( kind = 4 ) i 1732 integer ( kind = 4 ) ierror 1733 integer ( kind = 4 ) isgn 1734 integer ( kind = 4 ) istate 1735 integer ( kind = 4 ) ival 1736 integer ( kind = 4 ) length 1737 character ( len = * ) s 1738 1739 ! ************************************************************************* 1740 1741 ierror = 0 1742 istate = 0 1743 isgn = 1 1744 ival = 0 1745 1746 do i = 1, len_trim ( s ) 1747 1748 c = s(i:i) 1749 ! 1750 ! Haven't read anything. 1751 ! 1752 if ( istate == 0 ) then 1753 1754 if ( c == ' ' ) then 1755 1756 else if ( c == '-' ) then 1757 istate = 1 1758 isgn = -1 1759 else if ( c == '+' ) then 1760 istate = 1 1761 isgn = + 1 1762 else if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then 1763 istate = 2 1764 ival = ichar ( c ) - ichar ( '0' ) 1765 else 1766 ierror = 1 1767 return 1768 end if 1769 ! 1770 ! Have read the sign, expecting digits. 1771 ! 1772 else if ( istate == 1 ) then 1773 1774 if ( c == ' ' ) then 1775 1776 else if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then 1777 istate = 2 1778 ival = ichar ( c ) - ichar ( '0' ) 1779 else 1780 ierror = 1 1781 return 1782 end if 1783 ! 1784 ! Have read at least one digit, expecting more. 1785 ! 1786 else if ( istate == 2 ) then 1787 1788 if ( lle ( '0', c ) .and. lle ( c, '9' ) ) then 1789 ival = 10 * ival + ichar ( c ) - ichar ( '0' ) 1790 else 1791 ival = isgn * ival 1792 length = i - 1 1793 return 1794 end if 1795 1796 end if 1797 1798 end do 1799 ! 1800 ! If we read all the characters in the string, see if we're OK. 1801 ! 1802 if ( istate == 2 ) then 1803 ival = isgn * ival 1804 length = len_trim ( s ) 1805 else 1806 ierror = 1 1807 length = 0 1808 end if 1809 1810 return 1811 end subroutine s_to_i4
m_hamiltonian/s_to_r8 [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
s_to_r8
FUNCTION
.
INPUTS
OUTPUT
SOURCE
1827 subroutine s_to_r8 ( s, dval, ierror, length ) 1828 1829 !*****************************************************************************80 1830 ! 1831 !! S_TO_R8 reads an R8 from a string. 1832 ! 1833 ! Discussion: 1834 ! 1835 ! The routine will read as many characters as possible until it reaches 1836 ! the end of the string, or encounters a character which cannot be 1837 ! part of the number. 1838 ! 1839 ! Legal input is: 1840 ! 1841 ! 1 blanks, 1842 ! 2 '+' or '-' sign, 1843 ! 2.5 blanks 1844 ! 3 integer part, 1845 ! 4 decimal point, 1846 ! 5 fraction part, 1847 ! 6 'E' or 'e' or 'D' or 'd', exponent marker, 1848 ! 7 exponent sign, 1849 ! 8 exponent integer part, 1850 ! 9 exponent decimal point, 1851 ! 10 exponent fraction part, 1852 ! 11 blanks, 1853 ! 12 final comma or semicolon, 1854 ! 1855 ! with most quantities optional. 1856 ! 1857 ! Example: 1858 ! 1859 ! S DVAL 1860 ! 1861 ! '1' 1.0 1862 ! ' 1 ' 1.0 1863 ! '1A' 1.0 1864 ! '12,34,56' 12.0 1865 ! ' 34 7' 34.0 1866 ! '-1E2ABCD' -100.0 1867 ! '-1X2ABCD' -1.0 1868 ! ' 2E-1' 0.2 1869 ! '23.45' 23.45 1870 ! '-4.2E+2' -420.0 1871 ! '17d2' 1700.0 1872 ! '-14e-2' -0.14 1873 ! 'e2' 100.0 1874 ! '-12.73e-9.23' -12.73 * 10.0**(-9.23) 1875 ! 1876 ! Licensing: 1877 ! 1878 ! This code is distributed under the GNU LGPL license. 1879 ! 1880 ! Modified: 1881 ! 1882 ! 07 September 2004 1883 ! 1884 ! Author: 1885 ! 1886 ! John Burkardt 1887 ! 1888 ! Parameters: 1889 ! 1890 ! Input, character ( len = * ) S, the string containing the 1891 ! data to be read. Reading will begin at position 1 and 1892 ! terminate at the end of the string, or when no more 1893 ! characters can be read to form a legal real. Blanks, 1894 ! commas, or other nonnumeric data will, in particular, 1895 ! cause the conversion to halt. 1896 ! 1897 ! Output, real ( dp ) DVAL, the value read from the string. 1898 ! 1899 ! Output, integer ( kind = 4 ) IERROR, error flag. 1900 ! 0, no errors occurred. 1901 ! 1, 2, 6 or 7, the input number was garbled. The 1902 ! value of IERROR is the last type of input successfully 1903 ! read. For instance, 1 means initial blanks, 2 means 1904 ! a plus or minus sign, and so on. 1905 ! 1906 ! Output, integer ( kind = 4 ) LENGTH, the number of characters read 1907 ! to form the number, including any terminating 1908 ! characters such as a trailing comma or blanks. 1909 ! 1910 character c 1911 real ( dp ) dval 1912 integer ( kind = 4 ) ierror 1913 integer ( kind = 4 ) ihave 1914 integer ( kind = 4 ) isgn 1915 integer ( kind = 4 ) iterm 1916 integer ( kind = 4 ) jbot 1917 integer ( kind = 4 ) jsgn 1918 integer ( kind = 4 ) jtop 1919 integer ( kind = 4 ) length 1920 integer ( kind = 4 ) nchar 1921 integer ( kind = 4 ) ndig 1922 real ( dp ) rbot 1923 real ( dp ) rexp 1924 real ( dp ) rtop 1925 character ( len = * ) s 1926 1927 ! ************************************************************************* 1928 1929 nchar = len_trim ( s ) 1930 1931 ierror = 0 1932 dval = 0.0D+00 1933 length = -1 1934 isgn = 1 1935 rtop = 0 1936 rbot = 1 1937 jsgn = 1 1938 jtop = 0 1939 jbot = 1 1940 ihave = 1 1941 iterm = 0 1942 1943 do 1944 1945 length = length + 1 1946 1947 if ( nchar < length+1 ) then 1948 exit 1949 end if 1950 1951 c = s(length+1:length+1) 1952 ! 1953 ! Blank character. 1954 ! 1955 if ( c == ' ' ) then 1956 1957 if ( ihave == 2 ) then 1958 1959 else if ( ihave == 6 .or. ihave == 7 ) then 1960 iterm = 1 1961 else if ( 1 < ihave ) then 1962 ihave = 11 1963 end if 1964 ! 1965 ! Comma. 1966 ! 1967 else if ( c == ',' .or. c == ';' ) then 1968 1969 if ( ihave /= 1 ) then 1970 iterm = 1 1971 ihave = 12 1972 length = length + 1 1973 end if 1974 ! 1975 ! Minus sign. 1976 ! 1977 else if ( c == '-' ) then 1978 1979 if ( ihave == 1 ) then 1980 ihave = 2 1981 isgn = -1 1982 else if ( ihave == 6 ) then 1983 ihave = 7 1984 jsgn = -1 1985 else 1986 iterm = 1 1987 end if 1988 ! 1989 ! Plus sign. 1990 ! 1991 else if ( c == '+' ) then 1992 1993 if ( ihave == 1 ) then 1994 ihave = 2 1995 else if ( ihave == 6 ) then 1996 ihave = 7 1997 else 1998 iterm = 1 1999 end if 2000 ! 2001 ! Decimal point. 2002 ! 2003 else if ( c == '.' ) then 2004 2005 if ( ihave < 4 ) then 2006 ihave = 4 2007 else if ( 6 <= ihave .and. ihave <= 8 ) then 2008 ihave = 9 2009 else 2010 iterm = 1 2011 end if 2012 ! 2013 ! Scientific notation exponent marker. 2014 ! 2015 else if ( ch_eqi ( c, 'E' ) .or. ch_eqi ( c, 'D' ) ) then 2016 2017 if ( ihave < 6 ) then 2018 ihave = 6 2019 else 2020 iterm = 1 2021 end if 2022 ! 2023 ! Digit. 2024 ! 2025 else if ( ihave < 11 .and. lle ( '0', c ) .and. lle ( c, '9' ) ) then 2026 2027 if ( ihave <= 2 ) then 2028 ihave = 3 2029 else if ( ihave == 4 ) then 2030 ihave = 5 2031 else if ( ihave == 6 .or. ihave == 7 ) then 2032 ihave = 8 2033 else if ( ihave == 9 ) then 2034 ihave = 10 2035 end if 2036 2037 call ch_to_digit ( c, ndig ) 2038 2039 if ( ihave == 3 ) then 2040 rtop = 10.0D+00 * rtop + real ( ndig, dp ) 2041 else if ( ihave == 5 ) then 2042 rtop = 10.0D+00 * rtop + real ( ndig, dp ) 2043 rbot = 10.0D+00 * rbot 2044 else if ( ihave == 8 ) then 2045 jtop = 10 * jtop + ndig 2046 else if ( ihave == 10 ) then 2047 jtop = 10 * jtop + ndig 2048 jbot = 10 * jbot 2049 end if 2050 ! 2051 ! Anything else is regarded as a terminator. 2052 ! 2053 else 2054 iterm = 1 2055 end if 2056 ! 2057 ! If we haven't seen a terminator, and we haven't examined the 2058 ! entire string, go get the next character. 2059 ! 2060 if ( iterm == 1 ) then 2061 exit 2062 end if 2063 2064 end do 2065 ! 2066 ! If we haven't seen a terminator, and we have examined the 2067 ! entire string, then we're done, and LENGTH is equal to NCHAR. 2068 ! 2069 if ( iterm /= 1 .and. length+1 == nchar ) then 2070 length = nchar 2071 end if 2072 ! 2073 ! Number seems to have terminated. Have we got a legal number? 2074 ! Not if we terminated in states 1, 2, 6 or 7! 2075 ! 2076 if ( ihave == 1 .or. ihave == 2 .or. ihave == 6 .or. ihave == 7 ) then 2077 ierror = ihave 2078 write ( std_out, '(a)' ) ' ' 2079 write ( std_out, '(a)' ) 'S_TO_R8 - Serious error!' 2080 write ( std_out, '(a)' ) ' Illegal or nonnumeric input:' 2081 write ( std_out, '(a)' ) ' ' // trim ( s ) 2082 return 2083 end if 2084 ! 2085 ! Number seems OK. Form it. 2086 ! 2087 if ( jtop == 0 ) then 2088 rexp = 1.0D+00 2089 else 2090 if ( jbot == 1 ) then 2091 rexp = 10.0D+00 ** ( jsgn * jtop ) 2092 else 2093 rexp = 10.0D+00 ** ( real ( jsgn * jtop, dp ) & 2094 / real ( jbot, dp ) ) 2095 end if 2096 end if 2097 2098 dval = real ( isgn, dp ) * rexp * rtop / rbot 2099 2100 return 2101 end subroutine s_to_r8
m_hamiltonian/scqf [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
scqf
FUNCTION
.
INPUTS
OUTPUT
SOURCE
2117 subroutine scqf ( nt, t, mlt, wts, nwts, ndx, swts, st, gaussian_kind, alpha, beta, a, & 2118 b ) 2119 2120 !*****************************************************************************80 2121 ! 2122 !! SCQF scales a quadrature formula to a nonstandard interval. 2123 ! 2124 ! Discussion: 2125 ! 2126 ! The arrays WTS and SWTS may coincide. 2127 ! 2128 ! The arrays T and ST may coincide. 2129 ! 2130 ! Licensing: 2131 ! 2132 ! This code is distributed under the GNU LGPL license. 2133 ! 2134 ! Modified: 2135 ! 2136 ! 27 December 2009 2137 ! 2138 ! Author: 2139 ! 2140 ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. 2141 ! FORTRAN90 version by John Burkardt. 2142 ! 2143 ! Reference: 2144 ! 2145 ! Sylvan Elhay, Jaroslav Kautsky, 2146 ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of 2147 ! Interpolatory Quadrature, 2148 ! ACM Transactions on Mathematical Software, 2149 ! Volume 13, Number 4, December 1987, pages 399-415. 2150 ! 2151 ! Parameters: 2152 ! 2153 ! Input, integer ( kind = 4 ) NT, the number of knots. 2154 ! 2155 ! Input, real ( dp ) T(NT), the original knots. 2156 ! 2157 ! Input, integer ( kind = 4 ) MLT(NT), the multiplicity of the knots. 2158 ! 2159 ! Input, real ( dp ) WTS(NWTS), the weights. 2160 ! 2161 ! Input, integer ( kind = 4 ) NWTS, the number of weights. 2162 ! 2163 ! Input, integer ( kind = 4 ) NDX(NT), used to index the array WTS. 2164 ! For more details see the comments in CAWIQ. 2165 ! 2166 ! Output, real ( dp ) SWTS(NWTS), the scaled weights. 2167 ! 2168 ! Output, real ( dp ) ST(NT), the scaled knots. 2169 ! 2170 ! Input, integer ( kind = 4 ) KIND, the rule. 2171 ! 1, Legendre, (a,b) 1.0 2172 ! 2, Chebyshev Type 1, (a,b) ((b-x)*(x-a))^(-0.5) 2173 ! 3, Gegenbauer, (a,b) ((b-x)*(x-a))^alpha 2174 ! 4, Jacobi, (a,b) (b-x)^alpha*(x-a)^beta 2175 ! 5, Generalized Laguerre, (a,+oo) (x-a)^alpha*exp(-b*(x-a)) 2176 ! 6, Generalized Hermite, (-oo,+oo) |x-a|^alpha*exp(-b*(x-a)^2) 2177 ! 7, Exponential, (a,b) |x-(a+b)/2.0|^alpha 2178 ! 8, Rational, (a,+oo) (x-a)^alpha*(x+b)^beta 2179 ! 9, Chebyshev Type 2, (a,b) ((b-x)*(x-a))^(+0.5) 2180 ! 2181 ! Input, real ( dp ) ALPHA, the value of Alpha, if needed. 2182 ! 2183 ! Input, real ( dp ) BETA, the value of Beta, if needed. 2184 ! 2185 ! Input, real ( dp ) A, B, the interval endpoints. 2186 ! 2187 integer ( kind = 4 ) nt 2188 integer ( kind = 4 ) nwts 2189 2190 real ( dp ) a 2191 real ( dp ) al 2192 real ( dp ) alpha 2193 real ( dp ) b 2194 real ( dp ) be 2195 real ( dp ) beta 2196 integer ( kind = 4 ) i 2197 integer ( kind = 4 ) k 2198 integer ( kind = 4 ) gaussian_kind 2199 integer ( kind = 4 ) l 2200 integer ( kind = 4 ) mlt(nt) 2201 integer ( kind = 4 ) ndx(nt) 2202 real ( dp ) p 2203 real ( dp ) shft 2204 real ( dp ) slp 2205 real ( dp ) st(nt) 2206 real ( dp ) swts(nwts) 2207 real ( dp ) t(nt) 2208 real ( dp ) temp 2209 real ( dp ) tmp 2210 real ( dp ) wts(nwts) 2211 2212 ! ************************************************************************* 2213 2214 temp = epsilon ( temp ) 2215 2216 call parchk ( gaussian_kind, 1, alpha, beta ) 2217 2218 if ( gaussian_kind == 1 ) then 2219 2220 al = 0.0D+00 2221 be = 0.0D+00 2222 2223 if ( abs ( b - a ) <= temp ) then 2224 ABI_ERROR('SCQF - Fatal error! |B - A| too small.') 2225 end if 2226 2227 shft = ( a + b ) / 2.0D+00 2228 slp = ( b - a ) / 2.0D+00 2229 2230 else if ( gaussian_kind == 2 ) then 2231 2232 al = -0.5D+00 2233 be = -0.5D+00 2234 2235 if ( abs ( b - a ) <= temp ) then 2236 ABI_ERROR('SCQF - Fatal error! |B - A| too small.') 2237 end if 2238 2239 shft = ( a + b ) / 2.0D+00 2240 slp = ( b - a ) / 2.0D+00 2241 2242 else if ( gaussian_kind == 3 ) then 2243 2244 al = alpha 2245 be = alpha 2246 2247 if ( abs ( b - a ) <= temp ) then 2248 ABI_ERROR('SCQF - Fatal error! |B - A| too small.') 2249 end if 2250 2251 shft = ( a + b ) / 2.0D+00 2252 slp = ( b - a ) / 2.0D+00 2253 2254 else if ( gaussian_kind == 4 ) then 2255 2256 al = alpha 2257 be = beta 2258 2259 if ( abs ( b - a ) <= temp ) then 2260 ABI_ERROR("SCQF - Fatal error! |B - A| too small.") 2261 end if 2262 2263 shft = ( a + b ) / 2.0D+00 2264 slp = ( b - a ) / 2.0D+00 2265 2266 else if ( gaussian_kind == 5 ) then 2267 2268 if ( b <= 0.0D+00 ) then 2269 ABI_ERROR('SCQF - Fatal error! B <= 0') 2270 end if 2271 2272 shft = a 2273 slp = 1.0D+00 / b 2274 al = alpha 2275 be = 0.0D+00 2276 2277 else if ( gaussian_kind == 6 ) then 2278 2279 if ( b <= 0.0D+00 ) then 2280 ABI_ERROR('SCQF - Fatal error! B <= 0.') 2281 end if 2282 2283 shft = a 2284 slp = 1.0D+00 / sqrt ( b ) 2285 al = alpha 2286 be = 0.0D+00 2287 2288 else if ( gaussian_kind == 7 ) then 2289 2290 al = alpha 2291 be = 0.0D+00 2292 2293 if ( abs ( b - a ) <= temp ) then 2294 ABI_ERROR('|B - A| too small.') 2295 end if 2296 2297 shft = ( a + b ) / 2.0D+00 2298 slp = ( b - a ) / 2.0D+00 2299 2300 else if ( gaussian_kind == 8 ) then 2301 2302 if ( a + b <= 0.0D+00 ) then 2303 ABI_ERROR(' A + B <= 0.') 2304 end if 2305 2306 shft = a 2307 slp = a + b 2308 al = alpha 2309 be = beta 2310 2311 else if ( gaussian_kind == 9 ) then 2312 2313 al = 0.5D+00 2314 be = 0.5D+00 2315 2316 if ( abs ( b - a ) <= temp ) then 2317 ABI_ERROR('|B - A| too small.') 2318 end if 2319 2320 shft = ( a + b ) / 2.0D+00 2321 slp = ( b - a ) / 2.0D+00 2322 2323 end if 2324 2325 p = slp**( al + be + 1.0D+00 ) 2326 2327 do k = 1, nt 2328 2329 st(k) = shft + slp * t(k) 2330 l = abs ( ndx(k) ) 2331 2332 if ( l /= 0 ) then 2333 tmp = p 2334 do i = l, l + mlt(k) - 1 2335 swts(i) = wts(i) * tmp 2336 tmp = tmp * slp 2337 end do 2338 end if 2339 2340 end do 2341 2342 return 2343 end subroutine scqf
m_hamiltonian/sgqf [ Functions ]
[ Top ] [ m_hamiltonian ] [ Functions ]
NAME
sgqf
FUNCTION
.
INPUTS
OUTPUT
SOURCE
2360 subroutine sgqf ( nt, aj, bj, zemu, t, wts ) 2361 2362 !*****************************************************************************80 2363 ! 2364 !! SGQF computes knots and weights of a Gauss Quadrature formula. 2365 ! 2366 ! Discussion: 2367 ! 2368 ! This routine computes all the knots and weights of a Gauss quadrature 2369 ! formula with simple knots from the Jacobi matrix and the zero-th 2370 ! moment of the weight function, using the Golub-Welsch technique. 2371 ! 2372 ! Licensing: 2373 ! 2374 ! This code is distributed under the GNU LGPL license. 2375 ! 2376 ! Modified: 2377 ! 2378 ! 04 January 2010 2379 ! 2380 ! Author: 2381 ! 2382 ! Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky. 2383 ! FORTRAN90 version by John Burkardt. 2384 ! 2385 ! Reference: 2386 ! 2387 ! Sylvan Elhay, Jaroslav Kautsky, 2388 ! Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of 2389 ! Interpolatory Quadrature, 2390 ! ACM Transactions on Mathematical Software, 2391 ! Volume 13, Number 4, December 1987, pages 399-415. 2392 ! 2393 ! Parameters: 2394 ! 2395 ! Input, integer ( kind = 4 ) NT, the number of knots. 2396 ! 2397 ! Input, real ( dp ) AJ(NT), the diagonal of the Jacobi matrix. 2398 ! 2399 ! Input/output, real ( dp ) BJ(NT), the subdiagonal of the Jacobi 2400 ! matrix, in entries 1 through NT-1. On output, BJ has been overwritten. 2401 ! 2402 ! Input, real ( dp ) ZEMU, the zero-th moment of the weight function. 2403 ! 2404 ! Output, real ( dp ) T(NT), the knots. 2405 ! 2406 ! Output, real ( dp ) WTS(NT), the weights. 2407 ! 2408 integer ( kind = 4 ) nt 2409 2410 real ( dp ) aj(nt) 2411 real ( dp ) bj(nt) 2412 real ( dp ) t(nt) 2413 real ( dp ) wts(nt) 2414 real ( dp ) zemu 2415 2416 ! ************************************************************************* 2417 2418 ! 2419 ! Exit if the zero-th moment is not positive. 2420 ! 2421 if ( zemu <= 0.0D+00 ) then 2422 ABI_ERROR('SGQF - Fatal error! ZEMU <= 0.') 2423 end if 2424 ! 2425 ! Set up vectors for IMTQLX. 2426 ! 2427 t(1:nt) = aj(1:nt) 2428 2429 wts(1) = sqrt ( zemu ) 2430 wts(2:nt) = 0.0D+00 2431 ! 2432 ! Diagonalize the Jacobi matrix. 2433 ! 2434 call imtqlx ( nt, t, bj, wts ) 2435 2436 wts(1:nt) = wts(1:nt)**2 2437 2438 return 2439 end subroutine sgqf