GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis
111
|
Functions/Subroutines | |
subroutine, public | gl_expand_A_B (x, w, n, x_AB, w_AB, A, B) |
Takes the Gauss-Legendre rule for the interval [0,1] and expands it for the given interval [A,B]. More... | |
recursive real(kind=cfp) function, public | quad1d (f, A, B, eps, Qest) |
Adaptive 1D quadrature based on Gauss-Legendre rule. More... | |
recursive real(kind=cfp) function, public | quad2d (f, Ax, Bx, Ay, By, eps, Qest) |
Adaptive 2D quadrature on rectangle based on Gauss-Kronrod rule. The algorithm is based on that of Romanowski published in Int. J. Q. Chem. More... | |
real(kind=cfp) function, public | gl2d (f, Ax, Bx, Ay, By) |
2D Quadrature on rectangle using the Gauss-Kronrod rule of order 8. The meaning of the input variables is identical to quad2d input parameters. More... | |
subroutine, public | DQAGS (F, A, B, EPSABS, EPSREL, RESULT, ABSERR, NEVAL, IER, LIMIT, LENW, LAST, IWORK, WORK) |
***BEGIN PROLOGUE DQAGS ***PURPOSE The routine calculates an approximation result to a given Definite integral I = Integral of F over (A,B), Hopefully satisfying following claim for accuracy ABS(I-RESULT).LE.MAX(EPSABS,EPSREL*ABS(I)). ***LIBRARY SLATEC (QUADPACK) ***CATEGORY H2A1A1 ***TYPE real(kind=cfp) (QAGS-S, DQAGS-D) ***KEYWORDS AUTOMATIC INTEGRATOR, END POINT SINGULARITIES, EXTRAPOLATION, GENERAL-PURPOSE, GLOBALLY ADAPTIVE, QUADPACK, QUADRATURE ***AUTHOR Piessens, Robert Applied Mathematics and Programming Division K. U. Leuven de Doncker, Elise Applied Mathematics and Programming Division K. U. Leuven ***DESCRIPTION More... | |
subroutine, public | DQELG (N, EPSTAB, RESULT, ABSERR, RES3LA, NRES) |
***BEGIN PROLOGUE DQELG ***SUBSIDIARY ***PURPOSE The routine determines the limit of a given sequence of approximations, by means of the Epsilon algorithm of P.Wynn. An estimate of the absolute error is also given. The condensed Epsilon table is computed. Only those elements needed for the computation of the next diagonal are preserved. ***LIBRARY SLATEC ***TYPE real(kind=cfp) (QELG-S, DQELG-D) ***KEYWORDS CONVERGENCE ACCELERATION, EPSILON ALGORITHM, EXTRAPOLATION ***AUTHOR Piessens, Robert Applied Mathematics and Programming Division K. U. Leuven de Doncker, Elise Applied Mathematics and Programming Division K. U. Leuven ***DESCRIPTION More... | |
subroutine | DQK21 (F, A, B, RESULT, ABSERR, RESABS, RESASC) |
***BEGIN PROLOGUE DQK21 ***PURPOSE To compute I = Integral of F over (A,B), with error estimate J = Integral of ABS(F) over (A,B) ***LIBRARY SLATEC (QUADPACK) ***CATEGORY H2A1A2 ***TYPE real(kind=cfp) (QK21-S, DQK21-D) ***KEYWORDS 21-POINT GAUSS-KRONROD RULES, QUADPACK, QUADRATURE ***AUTHOR Piessens, Robert Applied Mathematics and Programming Division K. U. Leuven de Doncker, Elise Applied Mathematics and Programming Division K. U. Leuven ***DESCRIPTION More... | |
Variables | |
integer, parameter, public | n_7 = 7 |
Order of the Gauss-Legendre quadrature to which the x_7 and w_7 arrays correspond. More... | |
real(kind=cfp), dimension(2 *n_7+1), parameter, public | w_7 = (/0.015376620998058634177314196788602209_cfp,0.035183023744054062354633708225333669_cfp,0.05357961023358596750593477334293465_cfp,0.06978533896307715722390239725551416_cfp,0.08313460290849696677660043024060441_cfp,0.09308050000778110551340028093321141_cfp,0.09921574266355578822805916322191966_cfp,0.10128912096278063644031009998375966_cfp,0.09921574266355578822805916322191966_cfp,0.09308050000778110551340028093321141_cfp,0.08313460290849696677660043024060441_cfp,0.06978533896307715722390239725551416_cfp,0.05357961023358596750593477334293465_cfp,0.035183023744054062354633708225333669_cfp,0.015376620998058634177314196788602209_cfp/) |
Weights for the Gauss-Legendre quadrature of order 7 on interval [0,1]. More... | |
real(kind=cfp), dimension(2 *n_7+1), parameter, public | x_7 = (/0.0060037409897572857552171407066937094_cfp,0.031363303799647047846120526144895264_cfp,0.075896708294786391899675839612891574_cfp,0.13779113431991497629190697269303100_cfp,0.21451391369573057623138663137304468_cfp,0.30292432646121831505139631450947727_cfp,0.39940295300128273884968584830270190_cfp,0.50000000000000000000000000000000000_cfp,0.60059704699871726115031415169729810_cfp,0.69707567353878168494860368549052273_cfp,0.78548608630426942376861336862695532_cfp,0.86220886568008502370809302730696900_cfp,0.92410329170521360810032416038710843_cfp,0.96863669620035295215387947385510474_cfp,0.99399625901024271424478285929330629_cfp/) |
Abscissas for the Gauss-Legendre quadrature of order 7 on interval [0,1]. More... | |
integer, parameter, public | n_10 = 10 |
Order of the Gauss-Legendre quadrature to which the x_10 and w_10 arrays correspond. More... | |
real(kind=cfp), dimension(2 *n_10+1), parameter, public | x_10 = (/0.0031239146898052498698789820310295354_cfp,0.016386580716846852841688892546152419_cfp,0.039950332924799585604906433142515553_cfp,0.073318317708341358176374680706216165_cfp,0.11578001826216104569206107434688598_cfp,0.16643059790129384034701666500483042_cfp,0.22419058205639009647049060163784336_cfp,0.28782893989628060821316555572810597_cfp,0.35598934159879945169960374196769984_cfp,0.42721907291955245453148450883065683_cfp,0.50000000000000000000000000000000000_cfp,0.57278092708044754546851549116934317_cfp,0.64401065840120054830039625803230016_cfp,0.71217106010371939178683444427189403_cfp,0.77580941794360990352950939836215664_cfp,0.83356940209870615965298333499516958_cfp,0.88421998173783895430793892565311402_cfp,0.92668168229165864182362531929378384_cfp,0.96004966707520041439509356685748445_cfp,0.98361341928315314715831110745384758_cfp,0.99687608531019475013012101796897046_cfp/) |
Abscissas for the Gauss-Legendre quadrature of order 10 on interval [0,1]. More... | |
real(kind=cfp), dimension(2 *n_10+1), parameter, public | w_10 = (/0.008008614128887166662112308429235508_cfp,0.018476894885426246899975334149664833_cfp,0.028567212713428604141817913236223979_cfp,0.038050056814189651008525826650091590_cfp,0.046722211728016930776644870556966044_cfp,0.05439864958357418883173728903505282_cfp,0.06091570802686426709768358856286680_cfp,0.06613446931666873089052628724838780_cfp,0.06994369739553657736106671193379156_cfp,0.07226220199498502953191358327687627_cfp,0.07304056682484521359599257384168559_cfp,0.07226220199498502953191358327687627_cfp,0.06994369739553657736106671193379156_cfp,0.06613446931666873089052628724838780_cfp,0.06091570802686426709768358856286680_cfp,0.05439864958357418883173728903505282_cfp,0.046722211728016930776644870556966044_cfp,0.038050056814189651008525826650091590_cfp,0.028567212713428604141817913236223979_cfp,0.018476894885426246899975334149664833_cfp,0.008008614128887166662112308429235508_cfp/) |
Weights for the Gauss-Legendre quadrature of order 10 on interval [0,1]. More... | |
subroutine, public general_quadrature_gbl::DQAGS | ( | class(function_1d) | F, |
real(kind=cfp) | A, | ||
real(kind=cfp) | B, | ||
real(kind=cfp) | EPSABS, | ||
real(kind=cfp) | EPSREL, | ||
real(kind=cfp) | RESULT, | ||
real(kind=cfp) | ABSERR, | ||
integer | NEVAL, | ||
integer | IER, | ||
integer | LIMIT, | ||
integer | LENW, | ||
integer | LAST, | ||
integer, dimension(*) | IWORK, | ||
real(kind=cfp), dimension(*) | WORK | ||
) |
***BEGIN PROLOGUE DQAGS ***PURPOSE The routine calculates an approximation result to a given Definite integral I = Integral of F over (A,B), Hopefully satisfying following claim for accuracy ABS(I-RESULT).LE.MAX(EPSABS,EPSREL*ABS(I)). ***LIBRARY SLATEC (QUADPACK) ***CATEGORY H2A1A1 ***TYPE real(kind=cfp) (QAGS-S, DQAGS-D) ***KEYWORDS AUTOMATIC INTEGRATOR, END POINT SINGULARITIES, EXTRAPOLATION, GENERAL-PURPOSE, GLOBALLY ADAPTIVE, QUADPACK, QUADRATURE ***AUTHOR Piessens, Robert Applied Mathematics and Programming Division K. U. Leuven de Doncker, Elise Applied Mathematics and Programming Division K. U. Leuven ***DESCRIPTION
Computation of a definite integral Standard fortran subroutine Double precision version PARAMETERS ON ENTRY F - class(function_1d) Function whose method 'eval' defines the integrand Function F(X). A - Double precision Lower limit of integration B - Double precision Upper limit of integration EPSABS - Double precision Absolute accuracy requested EPSREL - Double precision Relative accuracy requested If EPSABS.LE.0 And EPSREL.LT.MAX(50*REL.MACH.ACC.,0.5D-28), The routine will end with IER = 6. ON RETURN RESULT - Double precision Approximation to the integral ABSERR - Double precision Estimate of the modulus of the absolute error, which should equal or exceed ABS(I-RESULT) NEVAL - Integer Number of integrand evaluations IER - Integer IER = 0 Normal and reliable termination of the routine. It is assumed that the requested accuracy has been achieved. IER.GT.0 Abnormal termination of the routine The estimates for integral and error are less reliable. It is assumed that the requested accuracy has not been achieved. ERROR MESSAGES IER = 1 Maximum number of subdivisions allowed has been achieved. One can allow more sub- divisions by increasing the value of LIMIT (and taking the according dimension adjustments into account. However, if this yields no improvement it is advised to analyze the integrand in order to determine the integration difficulties. If the position of a local difficulty can be determined (E.G. SINGULARITY, DISCONTINUITY WITHIN THE INTERVAL) one will probably gain from splitting up the interval at this point and calling the integrator on the subranges. If possible, an appropriate special-purpose integrator should be used, which is designed for handling the type of difficulty involved. = 2 The occurrence of roundoff error is detec- ted, which prevents the requested tolerance from being achieved. The error may be under-estimated. = 3 Extremely bad integrand behaviour occurs at some points of the integration interval. = 4 The algorithm does not converge. Roundoff error is detected in the Extrapolation table. It is presumed that the requested tolerance cannot be achieved, and that the returned result is the best which can be obtained. = 5 The integral is probably divergent, or slowly convergent. It must be noted that divergence can occur with any other value of IER. = 6 The input is invalid, because (EPSABS.LE.0 AND EPSREL.LT.MAX(50*REL.MACH.ACC.,0.5D-28) OR LIMIT.LT.1 OR LENW.LT.LIMIT*4. RESULT, ABSERR, NEVAL, LAST are set to zero. Except when LIMIT or LENW is invalid, IWORK(1), WORK(LIMIT*2+1) and WORK(LIMIT*3+1) are set to zero, WORK(1) is set to A and WORK(LIMIT+1) TO B. DIMENSIONING PARAMETERS LIMIT - Integer DIMENSIONING PARAMETER FOR IWORK LIMIT determines the maximum number of subintervals in the partition of the given integration interval (A,B), LIMIT.GE.1. IF LIMIT.LT.1, the routine will end with IER = 6. LENW - Integer DIMENSIONING PARAMETER FOR WORK LENW must be at least LIMIT*4. If LENW.LT.LIMIT*4, the routine will end with IER = 6. LAST - Integer On return, LAST equals the number of subintervals produced in the subdivision process, determines the number of significant elements actually in the WORK Arrays. WORK ARRAYS IWORK - Integer Vector of dimension at least LIMIT, the first K elements of which contain pointers to the error estimates over the subintervals such that WORK(LIMIT*3+IWORK(1)),... , WORK(LIMIT*3+IWORK(K)) form a decreasing sequence, with K = LAST IF LAST.LE.(LIMIT/2+2), and K = LIMIT+1-LAST otherwise WORK - Double precision Vector of dimension at least LENW on return WORK(1), ..., WORK(LAST) contain the left end-points of the subintervals in the partition of (A,B), WORK(LIMIT+1), ..., WORK(LIMIT+LAST) contain the right end-points, WORK(LIMIT*2+1), ..., WORK(LIMIT*2+LAST) contain the integral approximations over the subintervals, WORK(LIMIT*3+1), ..., WORK(LIMIT*3+LAST) contain the error estimates.
***REFERENCES (NONE) ***ROUTINES CALLED DQAGSE, XERMSG ***REVISION HISTORY (YYMMDD) 800101 DATE WRITTEN 890831 Modified array declarations. (WRB) 890831 REVISION DATE from Version 3.2 891214 Prologue converted to Version 4.0 format. (BAB) 900315 CALLs to XERROR changed to CALLs to XERMSG. (THJ) ***END PROLOGUE DQAGS
subroutine, public general_quadrature_gbl::DQELG | ( | integer | N, |
real(kind=cfp), dimension(max_epstab) | EPSTAB, | ||
real(kind=cfp) | RESULT, | ||
real(kind=cfp) | ABSERR, | ||
real(kind=cfp), dimension(3) | RES3LA, | ||
integer | NRES | ||
) |
***BEGIN PROLOGUE DQELG ***SUBSIDIARY ***PURPOSE The routine determines the limit of a given sequence of approximations, by means of the Epsilon algorithm of P.Wynn. An estimate of the absolute error is also given. The condensed Epsilon table is computed. Only those elements needed for the computation of the next diagonal are preserved. ***LIBRARY SLATEC ***TYPE real(kind=cfp) (QELG-S, DQELG-D) ***KEYWORDS CONVERGENCE ACCELERATION, EPSILON ALGORITHM, EXTRAPOLATION ***AUTHOR Piessens, Robert Applied Mathematics and Programming Division K. U. Leuven de Doncker, Elise Applied Mathematics and Programming Division K. U. Leuven ***DESCRIPTION
Epsilon algorithm Standard fortran subroutine Double precision version PARAMETERS N - Integer EPSTAB(N) contains the new element in the first column of the epsilon table. EPSTAB - Double precision Vector of dimension 52 containing the elements of the two lower diagonals of the triangular epsilon table. The elements are numbered starting at the right-hand corner of the triangle. RESULT - Double precision Resulting approximation to the integral ABSERR - Double precision Estimate of the absolute error computed from RESULT and the 3 previous results RES3LA - Double precision Vector of dimension 3 containing the last 3 results NRES - Integer Number of calls to the routine (should be zero at first call)
***SEE ALSO DQAGIE, DQAGOE, DQAGPE, DQAGSE ***ROUTINES CALLED F1MACH ***REVISION HISTORY (YYMMDD) 800101 DATE WRITTEN 890531 Changed all specific intrinsics to generic. (WRB) 890531 REVISION DATE from Version 3.2 891214 Prologue converted to Version 4.0 format. (BAB) 900328 Added TYPE section. (WRB) ***END PROLOGUE DQELG
LIST OF MAJOR VARIABLES ----------------------- E0 - THE 4 ELEMENTS ON WHICH THE COMPUTATION OF A NEW E1 ELEMENT IN THE EPSILON TABLE IS BASED E2 E3 E0 E3 E1 NEW E2 NEWELM - NUMBER OF ELEMENTS TO BE COMPUTED IN THE NEW DIAGONAL ERROR - ERROR = ABS(E1-E0)+ABS(E2-E1)+ABS(NEW-E2) RESULT - THE ELEMENT IN THE NEW DIAGONAL WITH LEAST VALUE OF ERROR MACHINE DEPENDENT CONSTANTS --------------------------- EPMACH IS THE LARGEST RELATIVE SPACING. OFLOW IS THE LARGEST POSITIVE MAGNITUDE. LIMEXP IS THE MAXIMUM NUMBER OF ELEMENTS THE EPSILON TABLE CAN CONTAIN. IF THIS NUMBER IS REACHED, THE UPPER DIAGONAL OF THE EPSILON TABLE IS DELETED.
subroutine general_quadrature_gbl::DQK21 | ( | class(function_1d) | F, |
real(kind=cfp) | A, | ||
real(kind=cfp) | B, | ||
real(kind=cfp) | RESULT, | ||
real(kind=cfp) | ABSERR, | ||
real(kind=cfp) | RESABS, | ||
real(kind=cfp) | RESASC | ||
) |
***BEGIN PROLOGUE DQK21 ***PURPOSE To compute I = Integral of F over (A,B), with error estimate J = Integral of ABS(F) over (A,B) ***LIBRARY SLATEC (QUADPACK) ***CATEGORY H2A1A2 ***TYPE real(kind=cfp) (QK21-S, DQK21-D) ***KEYWORDS 21-POINT GAUSS-KRONROD RULES, QUADPACK, QUADRATURE ***AUTHOR Piessens, Robert Applied Mathematics and Programming Division K. U. Leuven de Doncker, Elise Applied Mathematics and Programming Division K. U. Leuven ***DESCRIPTION
Integration rules Standard fortran subroutine Double precision version PARAMETERS ON ENTRY F - class(function_1d) Function whose method 'eval' defines the integrand Function F(X). A - Double precision Lower limit of integration B - Double precision Upper limit of integration ON RETURN RESULT - Double precision Approximation to the integral I RESULT is computed by applying the 21-POINT KRONROD RULE (RESK) obtained by optimal addition of abscissae to the 10-POINT GAUSS RULE (RESG). ABSERR - Double precision Estimate of the modulus of the absolute error, which should not exceed ABS(I-RESULT) RESABS - Double precision Approximation to the integral J RESASC - Double precision Approximation to the integral of ABS(F-I/(B-A)) over (A,B)
***REFERENCES (NONE) ***ROUTINES CALLED F1MACH ***REVISION HISTORY (YYMMDD) 800101 DATE WRITTEN 890531 Changed all specific intrinsics to generic. (WRB) 890531 REVISION DATE from Version 3.2 891214 Prologue converted to Version 4.0 format. (BAB) ***END PROLOGUE DQK21
THE ABSCISSAE AND WEIGHTS ARE GIVEN FOR THE INTERVAL (-1,1). BECAUSE OF SYMMETRY ONLY THE POSITIVE ABSCISSAE AND THEIR CORRESPONDING WEIGHTS ARE GIVEN. XGK - ABSCISSAE OF THE 21-POINT KRONROD RULE XGK(2), XGK(4), ... ABSCISSAE OF THE 10-POINT GAUSS RULE XGK(1), XGK(3), ... ABSCISSAE WHICH ARE OPTIMALLY ADDED TO THE 10-POINT GAUSS RULE WGK - WEIGHTS OF THE 21-POINT KRONROD RULE WG - WEIGHTS OF THE 10-POINT GAUSS RULE
GAUSS QUADRATURE WEIGHTS AND KRONROD QUADRATURE ABSCISSAE AND WEIGHTS AS EVALUATED WITH 80 DECIMAL DIGIT ARITHMETIC BY L. W. FULLERTON, BELL LABS, NOV. 1981.
LIST OF MAJOR VARIABLES ----------------------- CENTR - MID POINT OF THE INTERVAL HLGTH - HALF-LENGTH OF THE INTERVAL ABSC - ABSCISSA FVAL* - FUNCTION VALUE RESG - RESULT OF THE 10-POINT GAUSS FORMULA RESK - RESULT OF THE 21-POINT KRONROD FORMULA RESKH - APPROXIMATION TO THE MEAN VALUE OF F OVER (A,B), I.E. TO I/(B-A) MACHINE DEPENDENT CONSTANTS --------------------------- EPMACH IS THE LARGEST RELATIVE SPACING. UFLOW IS THE SMALLEST POSITIVE MAGNITUDE.
real(kind=cfp) function, public general_quadrature_gbl::gl2d | ( | class(function_2d) | f, |
real(kind=cfp), intent(in) | Ax, | ||
real(kind=cfp), intent(in) | Bx, | ||
real(kind=cfp), intent(in) | Ay, | ||
real(kind=cfp), intent(in) | By | ||
) |
2D Quadrature on rectangle using the Gauss-Kronrod rule of order 8. The meaning of the input variables is identical to quad2d input parameters.
subroutine, public general_quadrature_gbl::gl_expand_A_B | ( | real(kind=cfp), dimension(2*n+1), intent(in) | x, |
real(kind=cfp), dimension(2*n+1), intent(in) | w, | ||
integer, intent(in) | n, | ||
real(kind=cfp), dimension(2*n+1), intent(out) | x_AB, | ||
real(kind=cfp), dimension(2*n+1), intent(out) | w_AB, | ||
real(kind=cfp), intent(in) | A, | ||
real(kind=cfp), intent(in) | B | ||
) |
Takes the Gauss-Legendre rule for the interval [0,1] and expands it for the given interval [A,B].
recursive real(kind=cfp) function, public general_quadrature_gbl::quad1d | ( | class(function_1d) | f, |
real(kind=cfp), intent(in) | A, | ||
real(kind=cfp), intent(in) | B, | ||
real(kind=cfp), intent(in) | eps, | ||
real(kind=cfp), intent(in), optional | Qest | ||
) |
Adaptive 1D quadrature based on Gauss-Legendre rule.
[in] | f | The 1D function to be integrated. |
[in] | A | Interval start. |
[in] | B | Interval end. |
[in] | eps | Required relative precision for the integral. |
[in] | Qest | Optional: an estimate of the integral over the interval as obtained by a call to gl1d. |
recursive real(kind=cfp) function, public general_quadrature_gbl::quad2d | ( | class(function_2d) | f, |
real(kind=cfp), intent(in) | Ax, | ||
real(kind=cfp), intent(in) | Bx, | ||
real(kind=cfp), intent(in) | Ay, | ||
real(kind=cfp), intent(in) | By, | ||
real(kind=cfp), intent(in) | eps, | ||
real(kind=cfp), intent(in), optional | Qest | ||
) |
Adaptive 2D quadrature on rectangle based on Gauss-Kronrod rule. The algorithm is based on that of Romanowski published in Int. J. Q. Chem.
[in] | f | The 2D function to be integrated. |
[in] | Ax | Rectangle X-coordinate start. |
[in] | Bx | Rectangle X-coordinate end. |
[in] | Ay | Rectangle Y-coordinate start. |
[in] | By | Rectangle Y-coordinate end. |
[in] | eps | Required relative precision for the integral. |
[in] | Qest | Optional: an estimate of the integral over the specified rectangle as obtained by a call to gl2d. Note that the estimate must be obtained using gl2d for the algorithm to proceed correctly since it relies on the fixed-point G-K quadrature to calculate integrals on the sub-rectangles. |
integer, parameter, public general_quadrature_gbl::n_10 = 10 |
Order of the Gauss-Legendre quadrature to which the x_10 and w_10 arrays correspond.
integer, parameter, public general_quadrature_gbl::n_7 = 7 |
Order of the Gauss-Legendre quadrature to which the x_7 and w_7 arrays correspond.
real(kind=cfp), dimension(2*n_10+1), parameter, public general_quadrature_gbl::w_10 = (/0.008008614128887166662112308429235508_cfp,0.018476894885426246899975334149664833_cfp,0.028567212713428604141817913236223979_cfp,0.038050056814189651008525826650091590_cfp,0.046722211728016930776644870556966044_cfp,0.05439864958357418883173728903505282_cfp,0.06091570802686426709768358856286680_cfp,0.06613446931666873089052628724838780_cfp,0.06994369739553657736106671193379156_cfp,0.07226220199498502953191358327687627_cfp,0.07304056682484521359599257384168559_cfp,0.07226220199498502953191358327687627_cfp,0.06994369739553657736106671193379156_cfp,0.06613446931666873089052628724838780_cfp,0.06091570802686426709768358856286680_cfp,0.05439864958357418883173728903505282_cfp,0.046722211728016930776644870556966044_cfp,0.038050056814189651008525826650091590_cfp,0.028567212713428604141817913236223979_cfp,0.018476894885426246899975334149664833_cfp,0.008008614128887166662112308429235508_cfp/) |
Weights for the Gauss-Legendre quadrature of order 10 on interval [0,1].
real(kind=cfp), dimension(2*n_7+1), parameter, public general_quadrature_gbl::w_7 = (/0.015376620998058634177314196788602209_cfp,0.035183023744054062354633708225333669_cfp,0.05357961023358596750593477334293465_cfp,0.06978533896307715722390239725551416_cfp,0.08313460290849696677660043024060441_cfp,0.09308050000778110551340028093321141_cfp,0.09921574266355578822805916322191966_cfp,0.10128912096278063644031009998375966_cfp,0.09921574266355578822805916322191966_cfp,0.09308050000778110551340028093321141_cfp,0.08313460290849696677660043024060441_cfp,0.06978533896307715722390239725551416_cfp,0.05357961023358596750593477334293465_cfp,0.035183023744054062354633708225333669_cfp,0.015376620998058634177314196788602209_cfp/) |
Weights for the Gauss-Legendre quadrature of order 7 on interval [0,1].
real(kind=cfp), dimension(2*n_10+1), parameter, public general_quadrature_gbl::x_10 = (/0.0031239146898052498698789820310295354_cfp,0.016386580716846852841688892546152419_cfp,0.039950332924799585604906433142515553_cfp,0.073318317708341358176374680706216165_cfp,0.11578001826216104569206107434688598_cfp,0.16643059790129384034701666500483042_cfp,0.22419058205639009647049060163784336_cfp,0.28782893989628060821316555572810597_cfp,0.35598934159879945169960374196769984_cfp,0.42721907291955245453148450883065683_cfp,0.50000000000000000000000000000000000_cfp,0.57278092708044754546851549116934317_cfp,0.64401065840120054830039625803230016_cfp,0.71217106010371939178683444427189403_cfp,0.77580941794360990352950939836215664_cfp,0.83356940209870615965298333499516958_cfp,0.88421998173783895430793892565311402_cfp,0.92668168229165864182362531929378384_cfp,0.96004966707520041439509356685748445_cfp,0.98361341928315314715831110745384758_cfp,0.99687608531019475013012101796897046_cfp/) |
Abscissas for the Gauss-Legendre quadrature of order 10 on interval [0,1].
real(kind=cfp), dimension(2*n_7+1), parameter, public general_quadrature_gbl::x_7 = (/0.0060037409897572857552171407066937094_cfp,0.031363303799647047846120526144895264_cfp,0.075896708294786391899675839612891574_cfp,0.13779113431991497629190697269303100_cfp,0.21451391369573057623138663137304468_cfp,0.30292432646121831505139631450947727_cfp,0.39940295300128273884968584830270190_cfp,0.50000000000000000000000000000000000_cfp,0.60059704699871726115031415169729810_cfp,0.69707567353878168494860368549052273_cfp,0.78548608630426942376861336862695532_cfp,0.86220886568008502370809302730696900_cfp,0.92410329170521360810032416038710843_cfp,0.96863669620035295215387947385510474_cfp,0.99399625901024271424478285929330629_cfp/) |
Abscissas for the Gauss-Legendre quadrature of order 7 on interval [0,1].