ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_quadPack::setNodeWeightGK Interface Reference

Return the Kronrod \(2n + 1\) nodes and weights of the extension to the \(n\)-point Gauss-Legendre quadrature, as well as the \(n\)-point Gauss-Legendre weights. More...

Detailed Description

Return the Kronrod \(2n + 1\) nodes and weights of the extension to the \(n\)-point Gauss-Legendre quadrature, as well as the \(n\)-point Gauss-Legendre weights.

The abscissas and weights for both the Gauss and Kronrod rules are computed for integration over the interval \([-1,+1]\).
Also, only the nonnegative abscissas are computed because the quadrature formulae is symmetric with respect to the origin.
See the documentation of pm_quadPack for more details on the Gauss-Kronrod quadrature.

Parameters
[in]order: The input scalar integer of default kind IK containing the Number of Points in the Gauss quadrature rule.
(optional. If present, then the rest of the (output vector) arguments must be all allocatable.)
[out]nodeK: The output contiguous or allocatable vector of size (1 : n + 1) of type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128), containing the Gauss-Kronrod nodes (abscissas).
  1. If the input argument order is present, then Node must be an allocatable vector that will be reallocated to order+1 on output.
  2. If the input argument order is missing, then Node must be a contiguous vector of size order + 1, where order is the desired number of points in the Gauss rule.
[out]weightK: The output contiguous or allocatable vector of the same type, kind, and size as nodeK containing the Gauss-Kronrod optimal extension weights.
  1. If the input argument order is present, then weightK must be an allocatable vector that will be reallocated to order+1 on output.
  2. If the input argument order is missing, then weightK must be a contiguous vector of size order + 1, where order is the desired number of points in the Gauss rule.
[out]weightG: The output contiguous or allocatable vector of size (size(nodeK) / 2), of the same type and kind as nodeK containing the Gauss-Legendre weights.
By definition, the elements of weightG(:) correspond to the node elements nodeK(2::2).
  1. If the input argument order is present, then weightG must be an allocatable vector that will be reallocated to (order+1)/2 on output.
  2. If the input argument order is missing, then weightG must be a contiguous vector of size order + 1, where order is the desired number of points in the Gauss rule.


Possible calling interfaces

call setNodeWeightGK(nodeK, weightK, weightG)
Return the Kronrod nodes and weights of the extension to the -point Gauss-Legendre quadrature,...
This module contains classes and procedures for non-adaptive and adaptive global numerical quadrature...
Warning
The number of points \(n\) in the Gauss quadrature rule is determined from the size of the specified nodeK argument as n = size(nodeK) - 1.
As such, the sizes of weightK and weightG must be size(nodeK) - 1 and size(nodeK) / 2 respectively.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The procedures under discussion are impure.
See also
Robert Piessens, Maria Branders, 1974, A Note on the Optimal Addition of Abscissas to Quadrature Formulas of Gauss and Lobatto, Mathematics of Computation, 28, 125, 135-139.
Laurie, 1997, Calculation of gauss-kronrod quadrature rules.
Kronrod, 1965, Nodes and weights of quadrature formulas.
FORTRAN90 version by John Burkardt.


Example usage

1program example
2
3 use pm_kind, only: SK, IK, RKG => RKH ! testing with highest real precision available. all other real kinds are also supported.
4 use pm_quadpack, only: nodeK15, weightK15, weightG7
5 use pm_quadpack, only: nodeK21, weightK21, weightG10
6 use pm_quadpack, only: nodeK31, weightK31, weightG15
7 use pm_quadpack, only: nodeK41, weightK41, weightG20
8 use pm_quadpack, only: nodeK51, weightK51, weightG25
9 use pm_quadpack, only: nodeK61, weightK61, weightG30
10 use pm_quadpack, only: setNodeWeightGK
11 use pm_io, only: display_type
12
13 implicit none
14
15 integer(IK) , parameter :: MAX_NPG = 30_IK ! MAXimum Number of Points in the Gauss quadrature rule considered in this example.
16 real(RKG) :: nodeK(MAX_NPG+1), weightK(MAX_NPG+1), weightG(MAX_NPG+1)
17 real(RKG) , allocatable :: refNodeK(:), refWeightK(:), refWeightG(:)
18
19 type(display_type) :: disp
20 disp = display_type(file = "main.out.F90")
21
22 call disp%skip()
23 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
24 call disp%show("! Compute the 7-15 Gauss-Kronrod nodes and weights.")
25 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%skip()
27
28 refWeightG = real(weightG7, RKG)
29 refWeightK = real(weightK15, RKG)
30 refNodeK = real(nodeK15, RKG)
31 call setNodeWeight(npg = 7)
32
33 call disp%skip()
34 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
35 call disp%show("! Compute the 10-21 Gauss-Kronrod nodes and weights.")
36 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
37 call disp%skip()
38
39 refWeightG = real(weightG10, RKG)
40 refWeightK = real(weightK21, RKG)
41 refNodeK = real(nodeK21, RKG)
42 call setNodeWeight(npg = 10)
43
44 call disp%skip()
45 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
46 call disp%show("! Compute the 15-31 Gauss-Kronrod nodes and weights.")
47 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
48 call disp%skip()
49
50 refWeightG = real(weightG15, RKG)
51 refWeightK = real(weightK31, RKG)
52 refNodeK = real(nodeK31, RKG)
53 call setNodeWeight(npg = 15)
54
55 call disp%skip()
56 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
57 call disp%show("! Compute the 20-41 Gauss-Kronrod nodes and weights.")
58 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
59 call disp%skip()
60
61 refWeightG = real(weightG20, RKG)
62 refWeightK = real(weightK41, RKG)
63 refNodeK = real(nodeK41, RKG)
64 call setNodeWeight(npg = 20)
65
66 call disp%skip()
67 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
68 call disp%show("! Compute the 25-51 Gauss-Kronrod nodes and weights.")
69 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
70 call disp%skip()
71
72 refWeightG = real(weightG25, RKG)
73 refWeightK = real(weightK51, RKG)
74 refNodeK = real(nodeK51, RKG)
75 call setNodeWeight(npg = 25)
76
77 call disp%skip()
78 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
79 call disp%show("! Compute the 30-61 Gauss-Kronrod nodes and weights.")
80 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
81 call disp%skip()
82
83 refWeightG = real(weightG30, RKG)
84 refWeightK = real(weightK61, RKG)
85 refNodeK = real(nodeK61, RKG)
86 call setNodeWeight(npg = 30)
87
88contains
89
90 subroutine setNodeWeight(npg)
91 use pm_val2str, only: getStr
92 integer(IK) , intent(in) :: npg ! number of points in the Gauss quadrature rule.
93 integer(IK) :: npk ! number of extra points by the Kronrod extension rule ( = # points in the Gauss rule + 1).
94 integer(IK) :: i,j
95 call disp%skip()
96 call disp%show("npk = "//getStr(npg)//" + 1")
97 npk = npg + 1_IK
98 call disp%show("call setNodeWeightGK(nodeK(1:npk), weightK(1:npk), weightG(1:npk/2))")
99 call setNodeWeightGK(nodeK(1:npk), weightK(1:npk), weightG(1:npk/2))
100 call disp%show("[nodeK, weightK, weightG]")
101 do i = 1, npk
102 if (mod(i, 2_IK) == 0_IK) then
103 call disp%show( [nodeK(i), weightK(i), weightG(i/2)] )
104 else
105 call disp%show( [nodeK(i), weightK(i)] )
106 end if
107 end do
108 call disp%skip()
109 call disp%show("[nodeK - exact, weightK - exact, weightG - exact]")
110 do i = 1, npk
111 j = npk - i + 1
112 if (mod(i, 2_IK) == 0_IK) then
113 call disp%show( [nodeK(i) - refNodeK(j), weightK(i) - refWeightK(j), weightG(i/2) - refWeightG(size(refWeightG) - i/2 + 1)] )
114 else
115 call disp%show( [nodeK(i) - refNodeK(j), weightK(i) - refWeightK(j)] )
116 end if
117 end do
118 call disp%skip()
119 deallocate(refNodeK, refWeightG, refWeightK)
120 end subroutine
121
122end program example
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
Generate and return the conversion of the input value to an output Fortran string,...
Definition: pm_val2str.F90:167
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
real(RKB), dimension(15), parameter weightG30
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(21), parameter weightK41
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(16), parameter weightK31
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(21), parameter nodeK41
A vector of type real of the highest-precision kind available by the processor, containing the abscis...
real(RKB), dimension(8), parameter nodeK15
A vector of type real of the highest-precision kind available by the processor, containing the abscis...
real(RKB), dimension(11), parameter weightK21
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(31), parameter nodeK61
A vector of type real of the highest-precision kind available by the processor, containing the abscis...
real(RKB), dimension(16), parameter nodeK31
A vector of type real of the highest-precision kind available by the processor, containing the abscis...
real(RKB), dimension(5), parameter weightG10
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(13), parameter weightG25
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(11), parameter nodeK21
A vector of type real of the highest-precision kind available by the processor, containing the abscis...
real(RKB), dimension(8), parameter weightK15
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(4), parameter weightG7
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(26), parameter weightK51
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(10), parameter weightG20
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(26), parameter nodeK51
A vector of type real of the highest-precision kind available by the processor, containing the abscis...
real(RKB), dimension(31), parameter weightK61
A vector of type real of the highest-precision kind available by the processor, containing the weight...
real(RKB), dimension(8), parameter weightG15
A vector of type real of the highest-precision kind available by the processor, containing the weight...
This module contains the generic procedures for converting values of different types and kinds to For...
Definition: pm_val2str.F90:58
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3! Compute the 7-15 Gauss-Kronrod nodes and weights.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7npk = 7 + 1
8call setNodeWeightGK(nodeK(1:npk), weightK(1:npk), weightG(1:npk/2))
9[nodeK, weightK, weightG]
10+0.991455371120812639206854697526328534, +0.229353220105292249637320080589695566E-1
11+0.949107912342758524526189684047851219, +0.630920926299785532907006631892035206E-1, +0.129484966168869693270611432679081161
12+0.864864423359769072789712788640926167, +0.104790010322250183839876322541518041
13+0.741531185599394439863864773280788426, +0.140653259715525918745189590510237726, +0.279705391489276667901467771423779330
14+0.586087235467691130294144838258729596, +0.169004726639267902826583426598550208
15+0.405845151377397166906606412076961454, +0.190350578064785409913256402421013653, +0.381830050505118944950369775488975074
16+0.207784955007898467600689403773244930, +0.204432940075298892414161999234648970
17+0.00000000000000000000000000000000000, +0.209482141084727828012999174891714262, +0.417959183673469387755102040816326510
18
19[nodeK - exact, weightK - exact, weightG - exact]
20+0.00000000000000000000000000000000000, -0.361111864572606722447995864234673872E-34
21+0.00000000000000000000000000000000000, -0.770371977754894341222391177033970927E-33, -0.866668474974256133875190074163217293E-33
22+0.00000000000000000000000000000000000, +0.240741243048404481631997242823115915E-34
23+0.00000000000000000000000000000000000, -0.192592994438723585305597794258492732E-33, -0.240741243048404481631997242823115915E-33
24+0.00000000000000000000000000000000000, -0.722223729145213444895991728469347744E-34
25+0.00000000000000000000000000000000000, -0.240741243048404481631997242823115915E-34, -0.481482486096808963263994485646231830E-34
26+0.240741243048404481631997242823115915E-34, -0.120370621524202240815998621411557957E-33
27+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
28
29
30!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
31! Compute the 10-21 Gauss-Kronrod nodes and weights.
32!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
33
34
35npk = 10 + 1
36call setNodeWeightGK(nodeK(1:npk), weightK(1:npk), weightG(1:npk/2))
37[nodeK, weightK, weightG]
38+0.995657163025808080735527280689002850, +0.116946388673718742780643960621920455E-1
39+0.973906528517171720077964012084452104, +0.325581623079647274788189724593892905E-1, +0.666713443086881375935688098933313063E-1
40+0.930157491355708226001207180059508367, +0.547558965743519960313813002445801872E-1
41+0.865063366688984510732096688423493068, +0.750396748109199527670431409161898610E-1, +0.149451349150580593145776339657697151
42+0.780817726586416897063717578345042409, +0.931254545836976055350654650833663122E-1
43+0.679409568299024406234327365114873548, +0.109387158802297641899210590325804986, +0.219086362515982043995534934228163250
44+0.562757134668604683339000099272694155, +0.123491976262065851077958109831074157
45+0.433395394129247190799265943165784184, +0.134709217311473325928054001771706780, +0.269266719309996355091226921569469300
46+0.294392862701460198131126603103865508, +0.142775938577060080797094273138717120
47+0.148874338981631210884826001129719980, +0.147739104901338491374841515972067985, +0.295524224714752870173892994651338306
48+0.00000000000000000000000000000000000, +0.149445554002916905664936468389821279
49
50[nodeK - exact, weightK - exact, weightG - exact]
51+0.00000000000000000000000000000000000, -0.300926553810505602039996553528894894E-35
52+0.962964972193617926527988971292463659E-34, -0.469445423944388739182394623505076034E-33, -0.481482486096808963263994485646231830E-33
53+0.00000000000000000000000000000000000, +0.120370621524202240815998621411557957E-34
54+0.00000000000000000000000000000000000, -0.144444745829042688979198345693869549E-33, -0.192592994438723585305597794258492732E-33
55+0.00000000000000000000000000000000000, -0.361111864572606722447995864234673872E-34
56+0.00000000000000000000000000000000000, +0.240741243048404481631997242823115915E-34, +0.481482486096808963263994485646231830E-34
57+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
58+0.00000000000000000000000000000000000, -0.481482486096808963263994485646231830E-34, -0.481482486096808963263994485646231830E-34
59-0.481482486096808963263994485646231830E-34, +0.481482486096808963263994485646231830E-34
60+0.00000000000000000000000000000000000, -0.722223729145213444895991728469347744E-34, +0.00000000000000000000000000000000000
61+0.00000000000000000000000000000000000, +0.722223729145213444895991728469347744E-34
62
63
64!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
65! Compute the 15-31 Gauss-Kronrod nodes and weights.
66!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
67
68
69npk = 15 + 1
70call setNodeWeightGK(nodeK(1:npk), weightK(1:npk), weightG(1:npk/2))
71[nodeK, weightK, weightG]
72+0.998002298693397060285172840152271212, +0.537747987292334898779205143012766868E-2
73+0.987992518020485428489565718586612614, +0.150079473293161225383747630758065355E-1, +0.307532419961172683546283935772036556E-1
74+0.967739075679139134257347978784337170, +0.254608473267153201868740010196533869E-1
75+0.937273392400705904307758947710209522, +0.353463607913758462220379484783596456E-1, +0.703660474881081247092674164506669230E-1
76+0.897264532344081900882509656454495847, +0.445897513247648766082272993732796876E-1
77+0.848206583410427216200648320774216849, +0.534815246909280872653431472394303324E-1, +0.107159220467171935011869546685869334
78+0.790418501442465932967649294817947367, +0.620095678006706402851392309608029347E-1
79+0.724417731360170047416186054613938033, +0.698541213187282587095200770991473795E-1, +0.139570677926154314447804794511028232
80+0.650996741297416970533735895313274663, +0.768496807577203788944327774826590358E-1
81+0.570972172608538847537226737253910681, +0.830805028231330210382892472861039900E-1, +0.166269205816993933553200860481209062
82+0.485081863640239680693655740232350637, +0.885644430562117706472754436937743237E-1
83+0.394151347077563369897207370981045485, +0.931265981708253212254868727473457423E-1, +0.186161000015562211026800561866422892
84+0.299180007153168812166780024266388908, +0.966427269836236785051799076275893564E-1
85+0.201194093997434522300628303394596212, +0.991735987217919593323931734846031084E-1, +0.198431485327111576456118326443839338
86+0.101142066918717499027074231447392276, +0.100769845523875595044946662617569744
87+0.00000000000000000000000000000000000, +0.101330007014791549017374792767492326, +0.202578241925561272880620199967519315
88
89[nodeK - exact, weightK - exact, weightG - exact]
90+0.00000000000000000000000000000000000, +0.188079096131566001274997845955559308E-34
91+0.00000000000000000000000000000000000, -0.732756158528581140967391607842859066E-33, -0.761344181140579173161191280428104081E-33
92-0.962964972193617926527988971292463659E-34, +0.270833898429455041835996898176005404E-34
93+0.962964972193617926527988971292463659E-34, -0.403241582106077506733595381728719157E-33, -0.421297175334707842855995174940452851E-33
94+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
95+0.00000000000000000000000000000000000, +0.361111864572606722447995864234673872E-34, +0.361111864572606722447995864234673872E-34
96+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
97+0.00000000000000000000000000000000000, -0.962964972193617926527988971292463659E-34, -0.962964972193617926527988971292463659E-34
98+0.00000000000000000000000000000000000, +0.240741243048404481631997242823115915E-34
99+0.00000000000000000000000000000000000, +0.204630056591143809387197656399648528E-33, +0.240741243048404481631997242823115915E-33
100+0.00000000000000000000000000000000000, +0.240741243048404481631997242823115915E-34
101+0.00000000000000000000000000000000000, +0.240741243048404481631997242823115915E-34, +0.722223729145213444895991728469347744E-34
102-0.481482486096808963263994485646231830E-34, +0.240741243048404481631997242823115915E-34
103+0.00000000000000000000000000000000000, -0.240741243048404481631997242823115915E-34, +0.240741243048404481631997242823115915E-34
104-0.601853107621011204079993107057789787E-34, +0.240741243048404481631997242823115915E-34
105+0.00000000000000000000000000000000000, -0.216667118743564033468797518540804323E-33, +0.00000000000000000000000000000000000
106
107
108!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
109! Compute the 20-41 Gauss-Kronrod nodes and weights.
110!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
111
112
113npk = 20 + 1
114call setNodeWeightGK(nodeK(1:npk), weightK(1:npk), weightG(1:npk/2))
115[nodeK, weightK, weightG]
116+0.998859031588277663838315576545863027, +0.307358371852053150121829324603096036E-2
117+0.993128599185094924786122388471320311, +0.860026985564294219866178795010156739E-2, +0.176140071391521183118619623518520174E-1
118+0.981507877450250259193342994720216953, +0.146261692569712529837879603088683438E-1
119+0.963971927277913791267666131197277245, +0.203883734612665235980102314327544473E-1, +0.406014298003869413310399522749318474E-1
120+0.940822633831754753519982722212443384, +0.258821336049511588345050670961531257E-1
121+0.912234428251325905867752441203298155, +0.312873067770327989585431193238004526E-1, +0.626720483341090635695065351870413037E-1
122+0.878276811252281976077442995113078476, +0.366001697582007980305572407072109909E-1
123+0.839116971822218823394529061701520649, +0.416688733279736862637883059368949448E-1, +0.832767415767047487247581432220464148E-1
124+0.795041428837551198350638833272787951, +0.464348218674976747202318809261075123E-1
125+0.746331906460150792614305070355641587, +0.509445739237286919327076700503450663E-1, +0.101930119817240435036750135480350030
126+0.693237656334751384805490711845931522, +0.551951053482859947448323724197773487E-1
127+0.636053680726515025452836696226285966, +0.591114008806395723749672206485942313E-1, +0.118194531961518417312377377711382301
128+0.575140446819710315342946036586425083, +0.626532375547811680258701221742549883E-1
129+0.510867001950827098004364050955250976, +0.658345971336184221115635569693980091E-1, +0.131688638449176626898494499748163180
130+0.443593175238725103199992213492640121, +0.686486729285216193456234118853677711E-1
131+0.373706088715419560672548177024927261, +0.710544235534440683057903617232102466E-1, +0.142096109318382051329298325067165001
132+0.301627868114913004320555356858592298, +0.730306903327866674951894176589131234E-1
133+0.227785851141645078080496195368574615, +0.745828754004991889865814183624874938E-1, +0.149172986472603746787828737001969411
134+0.152605465240922675505220241022677516, +0.757044976845566746595427753766165649E-1
135+0.765265211334973337546404093988382039E-1, +0.763778676720807367055028350380610118E-1, +0.152753387130725850698084331955097591
136+0.00000000000000000000000000000000000, +0.766007119179996564450499015301016947E-1
137
138[nodeK - exact, weightK - exact, weightG - exact]
139+0.00000000000000000000000000000000000, -0.270833898429455041835996898176005404E-34
140+0.00000000000000000000000000000000000, -0.779399774369209509283591073639837774E-33, -0.800464633135944901426390832386860417E-33
141+0.00000000000000000000000000000000000, -0.120370621524202240815998621411557957E-34
142+0.00000000000000000000000000000000000, -0.258796836277034817754397036034849608E-33, -0.264815367353244929795196967105427506E-33
143+0.00000000000000000000000000000000000, -0.180555932286303361223997932117336936E-34
144+0.00000000000000000000000000000000000, -0.282870960581875265917596760317161200E-33, -0.300926553810505602039996553528894894E-33
145+0.00000000000000000000000000000000000, -0.180555932286303361223997932117336936E-34
146+0.00000000000000000000000000000000000, +0.204630056591143809387197656399648528E-33, +0.204630056591143809387197656399648528E-33
147+0.00000000000000000000000000000000000, -0.601853107621011204079993107057789787E-35
148+0.00000000000000000000000000000000000, +0.120370621524202240815998621411557957E-33, +0.156481807981462913060798207835025345E-33
149+0.00000000000000000000000000000000000, +0.180555932286303361223997932117336936E-34
150+0.00000000000000000000000000000000000, +0.120370621524202240815998621411557957E-34, +0.120370621524202240815998621411557957E-34
151-0.962964972193617926527988971292463659E-34, +0.120370621524202240815998621411557957E-34
152+0.00000000000000000000000000000000000, +0.601853107621011204079993107057789787E-34, +0.481482486096808963263994485646231830E-34
153+0.00000000000000000000000000000000000, -0.361111864572606722447995864234673872E-34
154+0.00000000000000000000000000000000000, +0.842594350669415685711990349880905702E-34, +0.722223729145213444895991728469347744E-34
155+0.481482486096808963263994485646231830E-34, +0.120370621524202240815998621411557957E-34
156+0.00000000000000000000000000000000000, -0.361111864572606722447995864234673872E-34, -0.240741243048404481631997242823115915E-34
157+0.00000000000000000000000000000000000, +0.120370621524202240815998621411557957E-34
158-0.120370621524202240815998621411557957E-34, +0.120370621524202240815998621411557957E-34, +0.00000000000000000000000000000000000
159+0.00000000000000000000000000000000000, -0.481482486096808963263994485646231830E-34
160
161
162!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
163! Compute the 25-51 Gauss-Kronrod nodes and weights.
164!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
165
166
167npk = 25 + 1
168call setNodeWeightGK(nodeK(1:npk), weightK(1:npk), weightG(1:npk/2))
169[nodeK, weightK, weightG]
170+0.999262104992609834193457486540340575, +0.198738389233031592650785188284342961E-2
171+0.995556969790498097908784946893901581, +0.556193213535671375804023690106697115E-2, +0.113937985010262879479029641132362507E-1
172+0.988035794534077247637331014577406252, +0.947397338617415160720771052365532276E-2
173+0.976663921459517511498315386479594089, +0.132362291955716748136564058469761715E-1, +0.263549866150321372619018152952990694E-1
174+0.961614986425842512418130033660167235, +0.168478177091282982315166675363363170E-1
175+0.942974571228974339414011169658470471, +0.204353711458828354565682922359383611E-1, +0.409391567013063126556234877116453113E-1
176+0.920747115281701561746346084546330642, +0.240099456069532162200924891648810722E-1
177+0.894991997878275368851042006782804903, +0.274753175878517378029484555178113998E-1, +0.549046959758351919259368915404736655E-1
178+0.865847065293275595448996969588340103, +0.307923001673874888911090202152285800E-1
179+0.833442628760834001421021108693569617, +0.340021302743293378367487952295509730E-1, +0.680383338123569172071871856567077122E-1
180+0.797873797998500059410410904994306606, +0.371162714834155435603306253676198537E-1
181+0.759259263037357630577282865204360970, +0.400838255040323820748392844670756497E-1, +0.801407003350010180132349596691112771E-1
182+0.717766406813084388186654079773297739, +0.428728450201700494768957924394951466E-1
183+0.673566368473468364485120633247622169, +0.455029130499217889098705847526603955E-1, +0.910282619829636498114972207028916491E-1
184+0.626810099010317412788122681624517907, +0.479825371388367139063922557569147330E-1
185+0.577662930241222967723689841612654042, +0.502776790807156719633252594334401640E-1, +0.100535949067050644202206890392685890
186+0.526325284334719182599623778158010128, +0.523628858064074758643667121378726962E-1
187+0.473002731445714960522182115009192053, +0.542511298885454901445433704598754928E-1, +0.108519624474263653116093957050116445
188+0.417885382193037748851814394594572502, +0.559508112204123173082406863827473103E-1
189+0.361172305809387837735821730127640676, +0.574371163615678328535826939395064512E-1, +0.114858259145711648339325545869555739
190+0.303089538931107830167478909980339298, +0.586896800223942079619741758567877471E-1
191+0.243866883720988432045190362797451579, +0.597203403241740599790992919325618536E-1, +0.119455763535784772228178126512901034
192+0.183718939421048892015969888759528396, +0.605394553760458629453602675175654328E-1
193+0.122864692610710396387359818808036808, +0.611285097170530483058590304162927646E-1, +0.122242442990310041688959518945851523
194+0.615444830056850788865463923667968797E-1, +0.614711898714253166615441319652640920E-1
195+0.00000000000000000000000000000000000, +0.615808180678329350787598242400646323E-1, +0.123176053726715451203902873079050139
196
197[nodeK - exact, weightK - exact, weightG - exact]
198+0.00000000000000000000000000000000000, +0.195602259976828641325997759793781681E-34
199+0.00000000000000000000000000000000000, +0.144896135659758447382258340524162891E-32, +0.147754937920958250601638307782687393E-32
200+0.00000000000000000000000000000000000, -0.150463276905252801019998276764447447E-35
201+0.00000000000000000000000000000000000, -0.662038418383112324487992417763568766E-34, -0.752316384526264005099991383822237234E-34
202+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
203-0.962964972193617926527988971292463659E-34, -0.613890169773431428161592969198945583E-33, -0.643982825154481988365592624551835072E-33
204+0.00000000000000000000000000000000000, -0.902779661431516806119989660586684681E-35
205-0.962964972193617926527988971292463659E-34, +0.321991412577240994182796312275917536E-33, +0.343056271343976386325596071022940179E-33
206+0.00000000000000000000000000000000000, -0.601853107621011204079993107057789787E-35
207+0.00000000000000000000000000000000000, -0.228704180895984257550397380681960119E-33, -0.252778305200824705713597104964271711E-33
208+0.00000000000000000000000000000000000, -0.240741243048404481631997242823115915E-34
209+0.00000000000000000000000000000000000, +0.601853107621011204079993107057789787E-35, -0.240741243048404481631997242823115915E-34
210+0.00000000000000000000000000000000000, -0.120370621524202240815998621411557957E-34
211+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
212+0.00000000000000000000000000000000000, -0.240741243048404481631997242823115915E-34
213+0.00000000000000000000000000000000000, +0.782409039907314565303991039175126723E-34, +0.601853107621011204079993107057789787E-34
214-0.962964972193617926527988971292463659E-34, -0.180555932286303361223997932117336936E-34
215+0.00000000000000000000000000000000000, -0.114352090447992128775198690340980060E-33, -0.180555932286303361223997932117336936E-33
216+0.00000000000000000000000000000000000, -0.361111864572606722447995864234673872E-34
217+0.00000000000000000000000000000000000, -0.180555932286303361223997932117336936E-34, -0.722223729145213444895991728469347744E-34
218-0.481482486096808963263994485646231830E-34, -0.180555932286303361223997932117336936E-34
219+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, -0.120370621524202240815998621411557957E-34
220-0.240741243048404481631997242823115915E-34, +0.601853107621011204079993107057789787E-35
221+0.00000000000000000000000000000000000, +0.541667796858910083671993796352010808E-34, +0.120370621524202240815998621411557957E-34
222+0.246759774124614593672797173893693813E-33, -0.842594350669415685711990349880905702E-34
223+0.00000000000000000000000000000000000, +0.782409039907314565303991039175126723E-34, +0.00000000000000000000000000000000000
224
225
226!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
227! Compute the 30-61 Gauss-Kronrod nodes and weights.
228!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
229
230
231npk = 30 + 1
232call setNodeWeightGK(nodeK(1:npk), weightK(1:npk), weightG(1:npk/2))
233[nodeK, weightK, weightG]
234+0.999484410050490637571325895705810795, +0.138901369867700762455159122675971330E-2
235+0.996893484074649540271630050918695337, +0.389046112709988405126720184451370613E-2, +0.796819249616660561546588347467177800E-2
236+0.991630996870404594858628366109485771, +0.663070391593129217331982636975016474E-2
237+0.983668123279747209970032581605662842, +0.927327965951776342844114689202402568E-2, +0.184664683110909591423021319120469131E-1
238+0.973116322501126268374693868423706846, +0.118230152534963417422328988532505948E-1
239+0.960021864968307512216871025581797646, +0.143697295070458048124514324435803472E-1, +0.287847078833233693497191796112923870E-1
240+0.944374444748559979415831324037439150, +0.169208891890532726275722894203221015E-1
241+0.926200047429274325879324277080474015, +0.194141411939423811734089510501285067E-1, +0.387991925696270495968019364463477550E-1
242+0.905573307699907798546522558925958353, +0.218280358216091922971674857383390032E-1
243+0.882560535792052681543116462530225587, +0.241911620780806013656863707252321658E-1, +0.484026728305940529029381404228076481E-1
244+0.857205233546061098958658510658943894, +0.265099548823331016106017093350754069E-1
245+0.829565762382768397442898119732501916, +0.287540487650412928439787853543341362E-1, +0.574931562176190664817216894020560505E-1
246+0.799727835821839083013668942322683280, +0.309072575623877624728842529430922737E-1
247+0.767777432104826194917977340974503127, +0.329814470574837260318141910168538974E-1, +0.659742298821804951281285151159623159E-1
248+0.733790062453226804726171131369527600, +0.349793380280600241374996707314678867E-1
249+0.697850494793315796932292388026640086, +0.368823646518212292239110656171359459E-1, +0.737559747377052062682438500221907194E-1
250+0.660061064126626961370053668149270727, +0.386789456247275929503486515322810469E-1
251+0.620526182989242861140477556431189339, +0.403745389515359591119952797524682789E-1, +0.807558952294202153546949384605298969E-1
252+0.579345235826361691756024932172540529, +0.419698102151642461471475412859697542E-1
253+0.536624148142019899264169793311072772, +0.434525397013560693168317281170732828E-1, +0.868997872010829798023875307151257365E-1
254+0.492480467861778574993693061207708789, +0.448148001331626631923555516167232373E-1
255+0.447033769538089176780609900322854008, +0.460592382710069881162717355593735825E-1, +0.921225222377861287176327070876187637E-1
256+0.400401254830394392535476211542660634, +0.471855465692991539452614781810994885E-1
257+0.352704725530878113471037207089373868, +0.481858617570871291407794922983045484E-1, +0.963687371746442596394686263518098187E-1
258+0.304073202273625077372677107199256562, +0.490554345550297788875281653672381693E-1
259+0.254636926167889846439805129817805114, +0.497956834270742063578115693799423231E-1, +0.995934205867952670627802821035694751E-1
260+0.204525116682309891438957671002024697, +0.504059214027823468408930856535850308E-1
261+0.153869913608583546963794672743255924, +0.508817958987496064922974730498046820E-1, +0.101762389748405504596428952168553985
262+0.102806937966737030147096751318000567, +0.512215478492587721706562826049441894E-1
263+0.514718425553176958330252131667225730E-1, +0.514261285374590259338628792157814533E-1, +0.102852652893558840341285636705415001
264+0.00000000000000000000000000000000000, +0.514947294294515675583404336470992344E-1
265
266[nodeK - exact, weightK - exact, weightG - exact]
267+0.00000000000000000000000000000000000, +0.135416949214727520917998449088002702E-34
268+0.962964972193617926527988971292463659E-34, -0.179728384263324470818387941595132475E-32, -0.184467977485839934050517887313212570E-32
269+0.00000000000000000000000000000000000, -0.376158192263132002549995691911118617E-35
270+0.00000000000000000000000000000000000, -0.334028474729661218264396174417073332E-33, -0.355093333496396610407195933164095974E-33
271+0.00000000000000000000000000000000000, +0.150463276905252801019998276764447447E-35
272+0.00000000000000000000000000000000000, +0.337037740267766274284796139952362281E-33, +0.343056271343976386325596071022940179E-33
273+0.00000000000000000000000000000000000, +0.902779661431516806119989660586684681E-35
274+0.00000000000000000000000000000000000, +0.511575141477859523467994140999121319E-34, +0.601853107621011204079993107057789787E-34
275+0.00000000000000000000000000000000000, +0.902779661431516806119989660586684681E-35
276+0.00000000000000000000000000000000000, +0.138426214752832576938398414623291651E-33, +0.132407683676622464897598483552713753E-33
277+0.00000000000000000000000000000000000, -0.601853107621011204079993107057789787E-35
278+0.00000000000000000000000000000000000, -0.752316384526264005099991383822237234E-34, -0.782409039907314565303991039175126723E-34
279+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
280+0.00000000000000000000000000000000000, -0.300926553810505602039996553528894894E-34, -0.481482486096808963263994485646231830E-34
281+0.00000000000000000000000000000000000, +0.120370621524202240815998621411557957E-34
282+0.00000000000000000000000000000000000, -0.240741243048404481631997242823115915E-34, -0.120370621524202240815998621411557957E-34
283+0.00000000000000000000000000000000000, -0.601853107621011204079993107057789787E-35
284+0.00000000000000000000000000000000000, +0.162500339057673025101598138905603243E-33, +0.168518870133883137142398069976181140E-33
285+0.00000000000000000000000000000000000, -0.601853107621011204079993107057789787E-35
286+0.00000000000000000000000000000000000, +0.240741243048404481631997242823115915E-34, +0.361111864572606722447995864234673872E-34
287+0.00000000000000000000000000000000000, -0.601853107621011204079993107057789787E-35
288+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
289+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
290+0.00000000000000000000000000000000000, -0.421297175334707842855995174940452851E-34, -0.481482486096808963263994485646231830E-34
291+0.00000000000000000000000000000000000, -0.601853107621011204079993107057789787E-35
292+0.00000000000000000000000000000000000, -0.601853107621011204079993107057789787E-35, +0.00000000000000000000000000000000000
293-0.240741243048404481631997242823115915E-34, +0.00000000000000000000000000000000000
294+0.00000000000000000000000000000000000, -0.120370621524202240815998621411557957E-34, -0.601853107621011204079993107057789787E-34
295-0.240741243048404481631997242823115915E-34, -0.180555932286303361223997932117336936E-34
296+0.00000000000000000000000000000000000, +0.192592994438723585305597794258492732E-33, -0.481482486096808963263994485646231830E-34
297+0.00000000000000000000000000000000000, -0.722223729145213444895991728469347744E-34
298
299
Test:
test_pm_quadPack


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, Oct 16, 2009, 11:14 AM, Michigan

Definition at line 2036 of file pm_quadPack.F90.


The documentation for this interface was generated from the following file: