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

Generate and return the approximate polynomial interpolation/extrapolation value of the input specified point x for the specified method. More...

Detailed Description

Generate and return the approximate polynomial interpolation/extrapolation value of the input specified point x for the specified method.

For polynomial interpolation/extrapolation, the computation relies on the Neville algorithm.
The extrapolation is done as if the out-of-bound query point is within the boundary (first or last) interpolation segment specified by the input (crdx, func) pairs of values.

Parameters
[in]method: The input scalar constant that can be,
  1. The scalar constant neimean or a scalar object of type neimean_type implying the use of the average of the func values of the two nearest neighbors of the input queryx smaller and larger than it as the output extrap.
  2. The scalar constant neinear or a scalar object of type neinear_type implying the use of the average of the func value of the neinear nearest neighbor of the input queryx as the output extrap.
    Note that the nearest neighbor in this case is measured by actual Euclidean distances of neighbors to the input queryx.
  3. The scalar constant neiprev or a scalar object of type neiprev_type implying the use of the func value of the largest abscissa in the input crdx smaller than the input queryx as the output extrap.
  4. The scalar constant neinext or a scalar object of type neinext_type implying the use of the func value of the smallest abscissa in the input crdx larger than the input queryx as the output extrap.
  5. The scalar constant piwilin or a scalar object of type piwilin_type implying the use of the linear interpolation/extrapolation of the func values of the two crdx points that bracket queryx as the output extrap.
    The linear interpolation/extrapolation implemented in this constructor is based on the Lagrange classical formula for linear interpolation/extrapolation.
    Suppose an input query point \(x\) falls between two nodes \(x_i\) and \(x_{i+1}\) with the corresponding function values \(y_i\) and \(y_{i+1}\) and we wish to estimate the corresponding interpolated/extrapolated value \(y(x)\), which can be computed as,

    \begin{equation*} y(x) = \frac {x - x_{i+1}} {x_i - x_{i+1}} y_i + \frac {x - x_{i}} {x_{i+1} - x_{i}} y_{i+1} ~. \end{equation*}

  6. The scalar constant monopol or a scalar object of type monopol_type implying the use of a single polynomial interpolation/extrapolation of highest degree size(crdx) - 1 possible to all pairs of (crdx, func) for computing the output extrap.
    The Neville algorithm is used to approximate the polynomial interpolation/extrapolation.
[in]crdx: The input contiguous vector of
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the set of abscissa in strictly ascending order through which the constructed polynomial must pass.
If needed, the input values for crdx and sortedY can be sorted in ascending order by calling either setSorted() or setSorted().
[in]func: The input contiguous vector of the same type, kind, and size as crdx, containing the set of function values corresponding to the input abscissa crdx through which the constructed polynomial must pass.
[in]queryx: The input scalar or vector of the same type and kind as crdx, containing the abscissa (x-value) of the queryx point.
Returns

extrap : The output object of the same type, kind, rank, and shape as queryx, containing the (approximate) interpolated/extrapolated y-value(s) corresponding to the queryx point(s).


Possible calling interfaces

extrap = getExtrap(method, crdx(1:nsam), func(1:nsam), queryx)
extrap(1:nque) = getExtrap(method, crdx(1:nsam), func(1:nsam), queryx(1:nque))
Generate and return the approximate polynomial interpolation/extrapolation value of the input specifi...
This module contains procedures and data types for interpolation of finite samples of data.
Warning
The condition size(crdx) == size(func) must hold for the corresponding input arguments.
The condition isAscendingAll(crdx) .or. same_type_as(method, monopol_type()) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
See also
getExtrap
setExtrap
getInterp
setInterp
pm_sampleQuan
pm_arraySort
pm_quadRomb


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
6 use pm_io, only: getErrTableWrite
7 use pm_io, only: display_type
8
9 implicit none
10
11 type(display_type) :: disp
12 integer(IK) :: itry, ntry = 5
13 integer(IK) :: dim, isam, ndim, nsam
14 disp = display_type(file = "main.out.F90")
15
16 call disp%skip()
17 call disp%show("!%%%%%%%%%%%%%%%%%")
18 call disp%show("!1D extrapolation.")
19 call disp%show("!%%%%%%%%%%%%%%%%%")
20 call disp%skip()
21
22 block
23 use pm_kind, only: RKG => RKS ! all processor kinds are supported.
24 use pm_mathConst, only: PI_RKH => PI
25 real(RKG), parameter :: PI = PI_RKH
26 real(RKG), allocatable :: crdx(:), func(:), queryx(:), extrap(:)
27 integer(IK) :: ncrd, nquery
28 call disp%skip()
29 call disp%show("ncrd = 9; nquery = 33")
30 ncrd = 9; nquery = 33
31 call disp%show("crdx = getLinSpace(0._RKG, 2 * PI, ncrd)")
32 crdx = getLinSpace(0._RKG, 2 * PI, ncrd)
33 call disp%show("crdx")
34 call disp%show( crdx )
35 call disp%show("func = sin(crdx)")
36 func = sin(crdx)
37 call disp%show("func")
38 call disp%show( func )
39 call disp%show("queryx = getLinSpace(-PI/2, 5 * PI/2, nquery)")
40 queryx = getLinSpace(-PI/2, 5 * PI/2, nquery)
41 call disp%show("queryx")
42 call disp%show( queryx )
43 call disp%skip()
44 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
45 call disp%show("!1D mean neighbors extrapolation.")
46 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
47 call disp%skip()
48 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.neimean.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'")
49 if (0 /= getErrTableWrite(SK_'getExtrap.neimean.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
50 call disp%show("extrap = getExtrap(neimean, crdx, func, queryx)")
51 extrap = getExtrap(neimean, crdx, func, queryx)
52 call disp%show("extrap")
53 call disp%show( extrap )
54 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.neimean.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'")
55 if (0 /= getErrTableWrite(SK_'getExtrap.neimean.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
56 call disp%skip()
57 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
58 call disp%show("!1D nearest neighbor extrapolation.")
59 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
60 call disp%skip()
61 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.neinear.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'")
62 if (0 /= getErrTableWrite(SK_'getExtrap.neinear.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
63 call disp%show("extrap = getExtrap(neinear, crdx, func, queryx)")
64 extrap = getExtrap(neinear, crdx, func, queryx)
65 call disp%show("extrap")
66 call disp%show( extrap )
67 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.neinear.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'")
68 if (0 /= getErrTableWrite(SK_'getExtrap.neinear.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
69 call disp%skip()
70 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
71 call disp%show("!1D next neighbor extrapolation.")
72 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
73 call disp%skip()
74 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.neinext.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'")
75 if (0 /= getErrTableWrite(SK_'getExtrap.neinext.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
76 call disp%show("extrap = getExtrap(neinext, crdx, func, queryx)")
77 extrap = getExtrap(neinext, crdx, func, queryx)
78 call disp%show("extrap")
79 call disp%show( extrap )
80 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.neinext.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'")
81 if (0 /= getErrTableWrite(SK_'getExtrap.neinext.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
82 call disp%skip()
83 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
84 call disp%show("!1D previous neighbor extrapolation.")
85 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
86 call disp%skip()
87 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.neiprev.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'")
88 if (0 /= getErrTableWrite(SK_'getExtrap.neiprev.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
89 call disp%show("extrap = getExtrap(neiprev, crdx, func, queryx)")
90 extrap = getExtrap(neiprev, crdx, func, queryx)
91 call disp%show("extrap")
92 call disp%show( extrap )
93 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.neiprev.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'")
94 if (0 /= getErrTableWrite(SK_'getExtrap.neiprev.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
95 call disp%skip()
96 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
97 call disp%show("!1D piecewise linear extrapolation.")
98 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
99 call disp%skip()
100 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.piwilin.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'")
101 if (0 /= getErrTableWrite(SK_'getExtrap.piwilin.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
102 call disp%show("extrap = getExtrap(piwilin, crdx, func, queryx)")
103 extrap = getExtrap(piwilin, crdx, func, queryx)
104 call disp%show("extrap")
105 call disp%show( extrap )
106 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.piwilin.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'")
107 if (0 /= getErrTableWrite(SK_'getExtrap.piwilin.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
108 call disp%skip()
109 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
110 call disp%show("!1D single polynomial extrapolation.")
111 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
112 call disp%skip()
113 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.monopol.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'")
114 if (0 /= getErrTableWrite(SK_'getExtrap.monopol.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
115 call disp%show("extrap = getExtrap(monopol, crdx, func, queryx)")
116 extrap = getExtrap(monopol, crdx, func, queryx)
117 call disp%show("extrap")
118 call disp%show( extrap )
119 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.monopol.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'")
120 if (0 /= getErrTableWrite(SK_'getExtrap.monopol.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
121 call disp%skip()
122 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
123 call disp%show("!1D single polynomial extrapolation showing Runge effect.")
124 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
125 call disp%skip()
126 block
127 use pm_distNorm, only: getNormLogPDF
128 call disp%show("ncrd = 9; nquery = 50 * ncrd")
129 ncrd = 9; nquery = 50 * ncrd
130 call disp%show("crdx = getLinSpace(-5., +5., ncrd)")
131 crdx = getLinSpace(-5., +5., ncrd)
132 call disp%show("crdx")
133 call disp%show( crdx )
134 call disp%show("func = exp(getNormLogPDF(crdx))")
135 func = exp(getNormLogPDF(crdx))
136 call disp%show("func")
137 call disp%show( func )
138 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.rungeEffect.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'")
139 if (0 /= getErrTableWrite(SK_'getExtrap.rungeEffect.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
140 call disp%show("queryx = getLinSpace(1.1 * minval(crdx, 1), 1.1 * maxval(crdx, 1), nquery)")
141 queryx = getLinSpace(1.1 * minval(crdx, 1), 1.1 * maxval(crdx, 1), nquery)
142 call disp%show("queryx")
143 call disp%show( queryx )
144 call disp%show("extrap = getExtrap(monopol, crdx, func, queryx)")
145 extrap = getExtrap(monopol, crdx, func, queryx)
146 call disp%show("extrap")
147 call disp%show( extrap )
148 call disp%show("if (0 /= getErrTableWrite(SK_'getExtrap.rungeEffect.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'")
149 if (0 /= getErrTableWrite(SK_'getExtrap.rungeEffect.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
150 call disp%skip()
151 end block
152 end block
153
154end program example
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Generate the natural logarithm of probability density function (PDF) of the univariate Normal distrib...
Generate and return the iostat code resulting from writing the input table of rank 1 or 2 to the spec...
Definition: pm_io.F90:5940
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
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
This module contains classes and procedures for computing various statistical quantities related to t...
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 LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
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 RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
This module contains relevant mathematical constants.
real(RKB), parameter PI
The scalar real constant of kind with highest available precision RKB representing the irrational num...
type(piwilin_type), parameter piwilin
type(neiprev_type), parameter neiprev
type(neinext_type), parameter neinext
type(monopol_type), parameter monopol
type(neinear_type), parameter neinear
type(neimean_type), parameter neimean
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!1D extrapolation.
4!%%%%%%%%%%%%%%%%%
5
6
7ncrd = 9; nquery = 33
8crdx = getLinSpace(0._RKG, 2 * PI, ncrd)
9crdx
10+0.00000000, +0.785398185, +1.57079637, +2.35619450, +3.14159274, +3.92699099, +4.71238899, +5.49778748, +6.28318548
11func = sin(crdx)
12func
13+0.00000000, +0.707106769, +1.00000000, +0.707106769, -0.874227766E-7, -0.707106888, -1.00000000, -0.707106531, +0.174845553E-6
14queryx = getLinSpace(-PI/2, 5 * PI/2, nquery)
15queryx
16-1.57079637, -1.27627206, -0.981747746, -0.687223434, -0.392699122, -0.981748104E-1, +0.196349502, +0.490873933, +0.785398126, +1.07992232, +1.37444675, +1.66897118, +1.96349537, +2.25801945, +2.55254412, +2.84706831, +3.14159250, +3.43611670, +3.73064089, +4.02516556, +4.31968975, +4.61421394, +4.90873861, +5.20326281, +5.49778700, +5.79231119, +6.08683538, +6.38136005, +6.67588472, +6.97040892, +7.26493311, +7.55945730, +7.85398197
17
18!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
19!1D mean neighbors extrapolation.
20!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
21
22if (0 /= getErrTableWrite(SK_'getExtrap.neimean.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
23extrap = getExtrap(neimean, crdx, func, queryx)
24extrap
25+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.353553385, +0.353553385, +0.353553385, +0.853553414, +0.853553414, +0.853553414, +0.853553414, +0.853553414, +0.353553355, +0.353553355, +0.353553355, -0.353553474, -0.353553474, -0.853553414, -0.853553414, -0.853553414, -0.853553295, -0.853553295, -0.853553295, -0.353553176, -0.353553176, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6
26if (0 /= getErrTableWrite(SK_'getExtrap.neimean.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
27
28!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
29!1D nearest neighbor extrapolation.
30!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
31
32if (0 /= getErrTableWrite(SK_'getExtrap.neinear.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
33extrap = getExtrap(neinear, crdx, func, queryx)
34extrap
35+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.707106769, +0.707106769, +0.707106769, +1.00000000, +1.00000000, +1.00000000, +0.707106769, +0.707106769, -0.874227766E-7, -0.874227766E-7, -0.874227766E-7, -0.707106888, -0.707106888, -0.707106888, -1.00000000, -1.00000000, -0.707106531, -0.707106531, -0.707106531, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6
36if (0 /= getErrTableWrite(SK_'getExtrap.neinear.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
37
38!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
39!1D next neighbor extrapolation.
40!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
41
42if (0 /= getErrTableWrite(SK_'getExtrap.neinext.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
43extrap = getExtrap(neinext, crdx, func, queryx)
44extrap
45+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.707106769, +0.707106769, +0.707106769, +1.00000000, +1.00000000, +0.707106769, +0.707106769, +0.707106769, -0.874227766E-7, -0.874227766E-7, -0.874227766E-7, -0.707106888, -0.707106888, -1.00000000, -1.00000000, -1.00000000, -0.707106531, -0.707106531, -0.707106531, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6
46if (0 /= getErrTableWrite(SK_'getExtrap.neinext.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
47
48!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
49!1D previous neighbor extrapolation.
50!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
51
52if (0 /= getErrTableWrite(SK_'getExtrap.neiprev.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
53extrap = getExtrap(neiprev, crdx, func, queryx)
54extrap
55+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.707106769, +0.707106769, +1.00000000, +1.00000000, +1.00000000, +0.707106769, +0.707106769, +0.707106769, -0.874227766E-7, -0.874227766E-7, -0.707106888, -0.707106888, -0.707106888, -1.00000000, -1.00000000, -1.00000000, -0.707106531, -0.707106531, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6, +0.174845553E-6
56if (0 /= getErrTableWrite(SK_'getExtrap.neiprev.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
57
58!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59!1D piecewise linear extrapolation.
60!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
61
62if (0 /= getErrTableWrite(SK_'getExtrap.piwilin.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
63extrap = getExtrap(piwilin, crdx, func, queryx)
64extrap
65-1.41421354, -1.14904845, -0.883883417, -0.618718445, -0.353553385, -0.883883759E-1, +0.176776648, +0.441941768, +0.707106709, +0.816941619, +0.926776648, +0.963388324, +0.853553414, +0.743718565, +0.530330002, +0.265165061, +0.127229356E-6, -0.265164793, -0.530329704, -0.743718445, -0.853553414, -0.963388264, -0.926776648, -0.816941679, -0.707106709, -0.441942036, -0.176777050, +0.883883536E-1, +0.353553772, +0.618718803, +0.883883774, +1.14904869, +1.41421402
66if (0 /= getErrTableWrite(SK_'getExtrap.piwilin.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
67
68!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
69!1D single polynomial extrapolation.
70!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
71
72if (0 /= getErrTableWrite(SK_'getExtrap.monopol.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
73extrap = getExtrap(monopol, crdx, func, queryx)
74extrap
75-1.78382504, -1.33006334, -0.990992546, -0.692502201, -0.398540974, -0.996969119E-1, +0.196270004, +0.472104400, +0.707106769, +0.881702781, +0.980677783, +0.995225310, +0.923959136, +0.773035169, +0.555531383, +0.290238500, +0.150949504E-6, -0.290238202, -0.555531144, -0.773034930, -0.923959076, -0.995225310, -0.980677783, -0.881702840, -0.707106888, -0.472104609, -0.196270570, +0.996969715E-1, +0.398542523, +0.692507029, +0.991004944, +1.33009934, +1.78385079
76if (0 /= getErrTableWrite(SK_'getExtrap.monopol.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
77
78!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
79!1D single polynomial extrapolation showing Runge effect.
80!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
81
82ncrd = 9; nquery = 50 * ncrd
83crdx = getLinSpace(-5., +5., ncrd)
84crdx
85-5.00000000, -3.75000000, -2.50000000, -1.25000000, +0.00000000, +1.25000000, +2.50000000, +3.75000000, +5.00000000
86func = exp(getNormLogPDF(crdx))
87func
88+0.148671938E-5, +0.352595642E-3, +0.175282992E-1, +0.182649091, +0.398942292, +0.182649091, +0.175282992E-1, +0.352595642E-3, +0.148671938E-5
89if (0 /= getErrTableWrite(SK_'getExtrap.rungeEffect.crd.txt', reshape([crdx, func], [ncrd, 2_IK]), header = SK_'crdx,func')) error stop 'extrapolation node outputting failed.'
90queryx = getLinSpace(1.1 * minval(crdx, 1), 1.1 * maxval(crdx, 1), nquery)
91queryx
92-5.50000000, -5.47550106, -5.45100212, -5.42650318, -5.40200424, -5.37750578, -5.35300684, -5.32850790, -5.30400896, -5.27951002, -5.25501108, -5.23051214, -5.20601320, -5.18151474, -5.15701580, -5.13251686, -5.10801792, -5.08351898, -5.05902004, -5.03452110, -5.01002216, -4.98552322, -4.96102428, -4.93652582, -4.91202688, -4.88752794, -4.86302900, -4.83853006, -4.81403112, -4.78953218, -4.76503325, -4.74053478, -4.71603584, -4.69153690, -4.66703796, -4.64253902, -4.61804008, -4.59354115, -4.56904221, -4.54454327, -4.52004433, -4.49554586, -4.47104692, -4.44654799, -4.42204905, -4.39755011, -4.37305117, -4.34855223, -4.32405376, -4.29955482, -4.27505589, -4.25055695, -4.22605801, -4.20155907, -4.17706013, -4.15256119, -4.12806225, -4.10356331, -4.07906437, -4.05456591, -4.03006697, -4.00556803, -3.98106909, -3.95657015, -3.93207121, -3.90757227, -3.88307357, -3.85857463, -3.83407593, -3.80957699, -3.78507805, -3.76057911, -3.73608017, -3.71158123, -3.68708253, -3.66258359, -3.63808465, -3.61358595, -3.58908701, -3.56458807, -3.54008913, -3.51559019, -3.49109149, -3.46659255, -3.44209361, -3.41759467, -3.39309573, -3.36859703, -3.34409809, -3.31959915, -3.29510021, -3.27060151, -3.24610257, -3.22160363, -3.19710469, -3.17260599, -3.14810705, -3.12360811, -3.09910917, -3.07461023, -3.05011153, -3.02561259, -3.00111365, -2.97661471, -2.95211601, -2.92761707, -2.90311813, -2.87861919, -2.85412025, -2.82962155, -2.80512261, -2.78062367, -2.75612473, -2.73162603, -2.70712709, -2.68262815, -2.65812922, -2.63363051, -2.60913157, -2.58463264, -2.56013370, -2.53563476, -2.51113605, -2.48663712, -2.46213818, -2.43763924, -2.41314054, -2.38864160, -2.36414266, -2.33964372, -2.31514478, -2.29064608, -2.26614714, -2.24164820, -2.21714926, -2.19265056, -2.16815162, -2.14365268, -2.11915374, -2.09465480, -2.07015610, -2.04565716, -2.02115822, -1.99665928, -1.97216058, -1.94766164, -1.92316270, -1.89866376, -1.87416506, -1.84966612, -1.82516718, -1.80066824, -1.77616930, -1.75167060, -1.72717166, -1.70267272, -1.67817378, -1.65367508, -1.62917614, -1.60467720, -1.58017826, -1.55567932, -1.53118062, -1.50668168, -1.48218298, -1.45768404, -1.43318510, -1.40868616, -1.38418722, -1.35968828, -1.33518934, -1.31069040, -1.28619146, -1.26169300, -1.23719406, -1.21269512, -1.18819618, -1.16369724, -1.13919830, -1.11469936, -1.09020042, -1.06570196, -1.04120302, -1.01670408, -0.992205143, -0.967706203, -0.943207264, -0.918708324, -0.894209385, -0.869710445, -0.845211983, -0.820713043, -0.796214104, -0.771715164, -0.747216225, -0.722717285, -0.698218346, -0.673719406, -0.649220467, -0.624722004, -0.600223064, -0.575724125, -0.551225185, -0.526726246, -0.502227306, -0.477728367, -0.453229427, -0.428730488, -0.404232025, -0.379733086, -0.355234146, -0.330735207, -0.306236267, -0.281737328, -0.257238388, -0.232739449, -0.208240509, -0.183742046, -0.159243107, -0.134744167, -0.110245228, -0.857462883E-1, -0.612473488E-1, -0.367484093E-1, -0.122494698E-1, +0.122494698E-1, +0.367479324E-1, +0.612468719E-1, +0.857458115E-1, +0.110244751, +0.134743690, +0.159242630, +0.183741570, +0.208240509, +0.232738972, +0.257237911, +0.281736851, +0.306235790, +0.330734730, +0.355233669, +0.379732609, +0.404231548, +0.428730488, +0.453228951, +0.477727890, +0.502226830, +0.526725769, +0.551224709, +0.575723648, +0.600222588, +0.624721527, +0.649220467, +0.673718929, +0.698217869, +0.722716808, +0.747215748, +0.771714687, +0.796213627, +0.820712566, +0.845211506, +0.869710445, +0.894208908, +0.918707848, +0.943206787, +0.967705727, +0.992204666, +1.01670361, +1.04120255, +1.06570148, +1.09020042, +1.11469889, +1.13919783, +1.16369677, +1.18819571, +1.21269464, +1.23719358, +1.26169252, +1.28619146, +1.31069040, +1.33518887, +1.35968781, +1.38418674, +1.40868568, +1.43318462, +1.45768356, +1.48218250, +1.50668144, +1.53117990, +1.55567884, +1.58017778, +1.60467672, +1.62917566, +1.65367460, +1.67817354, +1.70267248, +1.72717142, +1.75166988, +1.77616882, +1.80066776, +1.82516670, +1.84966564, +1.87416458, +1.89866352, +1.92316246, +1.94766140, +1.97215986, +1.99665880, +2.02115774, +2.04565668, +2.07015562, +2.09465456, +2.11915350, +2.14365244, +2.16815138, +2.19264984, +2.21714878, +2.24164772, +2.26614666, +2.29064560, +2.31514454, +2.33964348, +2.36414242, +2.38864136, +2.41313982, +2.43763876, +2.46213770, +2.48663664, +2.51113510, +2.53563404, +2.56013298, +2.58463192, +2.60913086, +2.63362980, +2.65812874, +2.68262768, +2.70712662, +2.73162556, +2.75612450, +2.78062344, +2.80512238, +2.82962132, +2.85412025, +2.87861919, +2.90311813, +2.92761707, +2.95211506, +2.97661400, +3.00111294, +3.02561188, +3.05011082, +3.07460976, +3.09910870, +3.12360764, +3.14810658, +3.17260551, +3.19710445, +3.22160339, +3.24610233, +3.27060127, +3.29510021, +3.31959915, +3.34409809, +3.36859608, +3.39309502, +3.41759396, +3.44209290, +3.46659184, +3.49109077, +3.51558971, +3.54008865, +3.56458759, +3.58908653, +3.61358547, +3.63808441, +3.66258335, +3.68708229, +3.71158123, +3.73608017, +3.76057911, +3.78507805, +3.80957603, +3.83407497, +3.85857391, +3.88307285, +3.90757179, +3.93207073, +3.95656967, +3.98106861, +4.00556755, +4.03006649, +4.05456543, +4.07906437, +4.10356331, +4.12806225, +4.15256119, +4.17706013, +4.20155907, +4.22605801, +4.25055599, +4.27505493, +4.29955387, +4.32405281, +4.34855175, +4.37305069, +4.39754963, +4.42204857, +4.44654751, +4.47104645, +4.49554539, +4.52004433, +4.54454327, +4.56904221, +4.59354115, +4.61804008, +4.64253902, +4.66703701, +4.69153595, +4.71603489, +4.74053383, +4.76503277, +4.78953171, +4.81403065, +4.83852959, +4.86302853, +4.88752747, +4.91202641, +4.93652534, +4.96102428, +4.98552322, +5.01002216, +5.03452110, +5.05902004, +5.08351898, +5.10801697, +5.13251591, +5.15701485, +5.18151379, +5.20601273, +5.23051167, +5.25501060, +5.27950954, +5.30400848, +5.32850742, +5.35300636, +5.37750530, +5.40200424, +5.42650318, +5.45100212, +5.47550106, +5.50000000
93extrap = getExtrap(monopol, crdx, func, queryx)
94extrap
95+1.84591424, +1.68530297, +1.53404069, +1.39173245, +1.25799561, +1.13246179, +1.01476717, +0.904567301, +0.801525354, +0.705316544, +0.615625620, +0.532148898, +0.454591900, +0.382672310, +0.316112578, +0.254649341, +0.198026896, +0.145998508, +0.983262807E-1, +0.547808707E-1, +0.151412832E-1, -0.208053477E-1, -0.532639138E-1, -0.824311003E-1, -0.108497709, -0.131645471, -0.152049571, -0.169878244, -0.185293078, -0.198448911, -0.209494323, -0.218571454, -0.225816920, -0.231360912, -0.235328153, -0.237837821, -0.239003673, -0.238934264, -0.237732947, -0.235498562, -0.232324779, -0.228300989, -0.223511547, -0.218037084, -0.211953789, -0.205333769, -0.198245347, -0.190752536, -0.182916731, -0.174794674, -0.166440085, -0.157903701, -0.149232596, -0.140471071, -0.131660253, -0.122838601, -0.114041790, -0.105302736, -0.966518819E-1, -0.881173611E-1, -0.797245502E-1, -0.714969933E-1, -0.634558946E-1, -0.556204244E-1, -0.480078310E-1, -0.406333320E-1, -0.335105062E-1, -0.266509764E-1, -0.200649593E-1, -0.137608796E-1, -0.774581823E-2, -0.202539330E-2, +0.339613250E-2, +0.851577707E-2, +0.133317299E-1, +0.178434551E-1, +0.220513791E-1, +0.259569362E-1, +0.295626335E-1, +0.328717791E-1, +0.358885676E-1, +0.386179537E-1, +0.410656184E-1, +0.432379991E-1, +0.451420993E-1, +0.467855111E-1, +0.481763557E-1, +0.493232496E-1, +0.502352864E-1, +0.509219170E-1, +0.513929985E-1, +0.516586453E-1, +0.517293103E-1, +0.516156703E-1, +0.513285846E-1, +0.508791469E-1, +0.502785482E-1, +0.495380759E-1, +0.486691482E-1, +0.476831794E-1, +0.465916470E-1, +0.454060547E-1, +0.441377275E-1, +0.427981243E-1, +0.413985252E-1, +0.399500579E-1, +0.384638533E-1, +0.369508006E-1, +0.354216918E-1, +0.338870808E-1, +0.323573612E-1, +0.308427252E-1, +0.293531120E-1, +0.278982297E-1, +0.264875069E-1, +0.251301378E-1, +0.238350090E-1, +0.226107500E-1, +0.214656480E-1, +0.204077084E-1, +0.194446128E-1, +0.185837168E-1, +0.178320482E-1, +0.171962772E-1, +0.166827533E-1, +0.162974615E-1, +0.160460528E-1, +0.159338061E-1, +0.159656405E-1, +0.161461290E-1, +0.164794717E-1, +0.169695076E-1, +0.176197160E-1, +0.184332132E-1, +0.194127560E-1, +0.205607153E-1, +0.218791477E-1, +0.233697221E-1, +0.250337496E-1, +0.268722028E-1, +0.288856979E-1, +0.310745575E-1, +0.334386714E-1, +0.359776653E-1, +0.386908054E-1, +0.415770747E-1, +0.446350649E-1, +0.478631519E-1, +0.512593016E-1, +0.548212640E-1, +0.585464910E-1, +0.624320805E-1, +0.664749369E-1, +0.706716031E-1, +0.750185251E-1, +0.795116946E-1, +0.841470212E-1, +0.889200345E-1, +0.938261747E-1, +0.988606289E-1, +0.104018286, +0.109293960, +0.114682131, +0.120177299, +0.125773579, +0.131465137, +0.137245759, +0.143109277, +0.149049297, +0.155059367, +0.161132768, +0.167262882, +0.173442885, +0.179665789, +0.185924858, +0.192212984, +0.198523000, +0.204847902, +0.211180508, +0.217513636, +0.223840088, +0.230152547, +0.236444131, +0.242707476, +0.248935461, +0.255121052, +0.261257023, +0.267336488, +0.273352474, +0.279298007, +0.285166204, +0.290950567, +0.296644330, +0.302241027, +0.307734221, +0.313117504, +0.318384796, +0.323530078, +0.328547269, +0.333430737, +0.338174939, +0.342774451, +0.347223759, +0.351518005, +0.355652094, +0.359621435, +0.363421291, +0.367047429, +0.370495617, +0.373761982, +0.376842737, +0.379734337, +0.382433534, +0.384937078, +0.387242287, +0.389346391, +0.391247064, +0.392942011, +0.394429386, +0.395707428, +0.396774769, +0.397630036, +0.398272455, +0.398701042, +0.398915499, +0.398915440, +0.398701042, +0.398272395, +0.397630036, +0.396774739, +0.395707428, +0.394429386, +0.392942011, +0.391247034, +0.389346421, +0.387242317, +0.384937137, +0.382433563, +0.379734427, +0.376842797, +0.373762071, +0.370495677, +0.367047459, +0.363421410, +0.359621525, +0.355652213, +0.351518095, +0.347223908, +0.342774481, +0.338175088, +0.333430886, +0.328547299, +0.323530108, +0.318384886, +0.313117594, +0.307734281, +0.302241147, +0.296644449, +0.290950626, +0.285166293, +0.279297978, +0.273352593, +0.267336607, +0.261257172, +0.255121171, +0.248935610, +0.242707640, +0.236444280, +0.230152696, +0.223840073, +0.217513740, +0.211180627, +0.204848036, +0.198523134, +0.192213073, +0.185924992, +0.179665923, +0.173442915, +0.167262882, +0.161132872, +0.155059442, +0.149049431, +0.143109396, +0.137245879, +0.131465226, +0.125773668, +0.120177366, +0.114682302, +0.109294072, +0.104018405, +0.988607034E-1, +0.938262716E-1, +0.889201090E-1, +0.841470733E-1, +0.795117617E-1, +0.750185698E-1, +0.706717297E-1, +0.664750114E-1, +0.624321550E-1, +0.585465655E-1, +0.548213422E-1, +0.512593761E-1, +0.478631854E-1, +0.446351022E-1, +0.415770933E-1, +0.386908874E-1, +0.359777212E-1, +0.334387198E-1, +0.310745947E-1, +0.288857389E-1, +0.268722363E-1, +0.250337645E-1, +0.233697332E-1, +0.218791589E-1, +0.205607507E-1, +0.194127727E-1, +0.184332393E-1, +0.176197328E-1, +0.169695131E-1, +0.164794754E-1, +0.161461271E-1, +0.159656331E-1, +0.159338024E-1, +0.160460435E-1, +0.162974577E-1, +0.166827403E-1, +0.171962641E-1, +0.178320222E-1, +0.185836945E-1, +0.194445904E-1, +0.204076767E-1, +0.214656107E-1, +0.226107091E-1, +0.238349903E-1, +0.251301117E-1, +0.264874734E-1, +0.278981999E-1, +0.293531083E-1, +0.308427252E-1, +0.323573612E-1, +0.338870659E-1, +0.354216844E-1, +0.369508043E-1, +0.384638608E-1, +0.399501026E-1, +0.413984731E-1, +0.427981019E-1, +0.441377461E-1, +0.454059988E-1, +0.465917066E-1, +0.476831943E-1, +0.486690998E-1, +0.495380983E-1, +0.502785295E-1, +0.508791432E-1, +0.513285920E-1, +0.516156778E-1, +0.517292917E-1, +0.516586602E-1, +0.513930023E-1, +0.509219021E-1, +0.502352864E-1, +0.493232906E-1, +0.481763929E-1, +0.467855334E-1, +0.451421440E-1, +0.432380587E-1, +0.410657004E-1, +0.386179984E-1, +0.358886197E-1, +0.328718424E-1, +0.295627043E-1, +0.259570107E-1, +0.220514089E-1, +0.178434961E-1, +0.133317783E-1, +0.851577334E-2, +0.339613203E-2, -0.202539284E-2, -0.774581730E-2, -0.137606375E-1, -0.200647190E-1, -0.266507715E-1, -0.335103124E-1, -0.406331867E-1, -0.480076671E-1, -0.556202531E-1, -0.634557307E-1, -0.714968443E-1, -0.797243416E-1, -0.881172419E-1, -0.966519117E-1, -0.105302691, -0.114041805, -0.122838676, -0.131660223, -0.140471101, -0.149232566, -0.157903314, -0.166439772, -0.174794197, -0.182916522, -0.190752506, -0.198245049, -0.205333620, -0.211953610, -0.218036979, -0.223511428, -0.228300855, -0.232324779, -0.235498562, -0.237732947, -0.238934264, -0.239003673, -0.237837821, -0.235328302, -0.231361106, -0.225817204, -0.218571752, -0.209494501, -0.198449165, -0.185293362, -0.169878602, -0.152049929, -0.131645873, -0.108498201, -0.824316442E-1, -0.532639138E-1, -0.208053477E-1, +0.151412832E-1, +0.547808707E-1, +0.983262807E-1, +0.145998508, +0.198024780, +0.254647046, +0.316110075, +0.382669598, +0.454590619, +0.532147288, +0.615623891, +0.705314636, +0.801523566, +0.904565215, +1.01476502, +1.13245940, +1.25799561, +1.39173245, +1.53404069, +1.68530297, +1.84591424
96if (0 /= getErrTableWrite(SK_'getExtrap.rungeEffect.extrap.txt', reshape([queryx, extrap], [nquery, 2_IK]), header = SK_'queryx,extrap')) error stop 'extrapolation outputting failed.'
97
98

Postprocessing of the example output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6import glob
7import sys
8
9linewidth = 2
10fontsize = 17
11
12for kind in ["neimean", "neinear", "neinext", "neiprev", "piwilin", "monopol", "rungeEffect"]:
13
14 crd = pd.read_csv(glob.glob("*."+kind+".crd.txt")[0], delimiter = ",")
15 pattern = "*."+kind+".extrap.txt"
16 fileList = glob.glob(pattern)
17
18 for file in fileList:
19
20 df = pd.read_csv(file, delimiter = ",")
21
22 # start with a square Figure
23 fig = plt.figure(figsize = (8, 6))
24 ax = plt.subplot(1,1,1)
25
26 ax.scatter ( df.values[:, 0]
27 , df.values[:, 1]
28 , zorder = 1000
29 , c = "black"
30 , s = 8
31 )
32 ax.scatter ( crd.values[:, 0]
33 , crd.values[:,1]
34 , zorder = 1000
35 , c = "red"
36 , s = 20
37 )
38
39 plt.minorticks_on()
40 ax.set_xlabel("X", fontsize = 17)
41 ax.set_ylabel("Y", fontsize = 17)
42 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
43 ax.tick_params(axis = "x", which = "minor")
44 ax.tick_params(axis = "y", which = "minor")
45 ax.legend([file.split(".")[-3], "nodes"], fontsize = fontsize)
46
47 plt.tight_layout()
48 plt.savefig(file.replace(".txt",".png"))

Visualization of the example output
Test:
test_pm_polation


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, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin

Definition at line 2540 of file pm_polation.F90.


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