Generate and return .true.
if the parameters of a least-squares fit to the histogram representing a Cyclic-Geometric-distributed sample can be successfully inferred, otherwise, return .false.
.
See the documentation of pm_distGeomCyclic for more details on the Cyclic Geometric distribution.
Given a set of Cyclic Bernoulli trial success steps \(\ms{stepSuccess}\) and the corresponding number of successes \(\ms{freqSuccess}\) at these steps, with the supplied Cyclic Bernoulli trial \(\ms{period}\), the procedures of this generic interface return the best-fit parameters \((\ms{probSuccess}, \ms{normFac})\) of the following curve,
\begin{eqnarray}
\large
\ms{freqSuccess}
&=& \ms{normFac} \times \pi_{\mathcal{CG}} (X = \ms{stepSuccess} ~|~ \ms{probSuccess}, \ms{period}) ~, \nonumber \\
&=& \ms{normFac} \times \frac{\ms{probSuccess} (1 - \ms{probSuccess})^{\ms{stepSuccess} - 1}}{1 - (1 - \ms{probSuccess})^{\ms{period}}} ~,
\end{eqnarray}
where the probability of success on each trial is \(\ms{probSuccess}\), and \(\pi_{\mathcal{CG}}\) is the probability mass function of the \(\ms{stepSuccess}\)th trial being the first success in a cyclic trial set with the cycle \(\ms{period}\).
- Note
- This fitting appears frequently in Fork-Join parallel MCMC simulations of the ParaMonte library.
See Amir Shahmoradi, Fatemeh Bagheri (2020). ParaDRAM: A Cross-Language Toolbox for Parallel High-Performance Delayed-Rejection Adaptive Metropolis Markov Chain Monte Carlo Simulations. for details and discussion.
- Parameters
-
[in] | stepSuccess | : The input contiguous positive vector of
-
type
real of the same kind as that of probSuccess ,
-
type
integer of default kind IK,
containing the steps ( \(x\)) at which contribution frequencies represented by input vector freqSuccess have occurred.
|
[in] | freqSuccess | : The input contiguous positive vector of the same type, kind, and size as the input stepSuccess , containing the contribution frequencies at the corresponding stepSuccess in the histogram.
|
[in] | period | : The input positive scalar of type integer of default kind IK, containing the period of the Cyclic Geometric distribution.
The period represents the maximum number of Bernoulli trials in the experiment.
Note that the condition maxval(stepSuccess) <= period must hold and the sizes of stepSuccess and freqSuccess vectors must be smaller than the specified period .
|
[out] | probSuccess | : The output positive scalar of
-
type
real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
representing the best-fit probability-of-success parameter of the Cyclic Geometric distribution inferred from the fitting.
|
[in] | normFac | : The output positive scalar of the same type and kind as the output probSuccess , containing the best-fit normalization constant of the histogram represented by the input stepSuccess and freqSuccess`.
|
- Returns
failed
: The output scalar of type logical
of default kind LK that is .true.
if and only if the best-fit output parameters are successfully inferred.
The algorithm can fail only if the optimizer fails to converge to set of best-fit parameters.
Possible calling interfaces ⛓
failed
= isFailedGeomCyclicFit(stepSuccess(
1:period), freqSuccess(
1:period), period, probSuccess, normFac)
Generate and return .true. if the parameters of a least-squares fit to the histogram representing a C...
This module contains classes and procedures for computing various statistical quantities related to t...
- Warning
- The condition
all(0 < stepSuccess)
must hold for the corresponding input arguments.
The condition 1 < size(freqSuccess))
must hold for the corresponding input arguments.
The condition size(stepSuccess) == size(freqSuccess)
must hold for the corresponding input arguments.
The condition all(stepSuccess <= size(stepSuccess))
must hold for the corresponding input arguments.
The condition all(stepSuccess <= period)
must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1
.
- See also
- setGeomCyclicRand
setGeomCyclicLogPMF
setGeomCyclicLogCDF
Example usage ⛓
14 integer(IK) :: period, i
15 type(display_type) :: disp
16 integer(IK),
allocatable :: rand(:)
17 integer(IK),
allocatable :: freqSuccess(:)
18 integer(IK),
allocatable :: stepSuccess(:)
23 real(RKG) :: probSuccess, probSuccessFit, normFacFit
24 real(RKG),
allocatable :: x(:), y(:)
26 call disp%show(
"probSuccess = .1; period = 20")
27 probSuccess
= .
1; period
= 20
28 call disp%show(
"rand = getGeomCyclicRand(probSuccess, [(period, i = 1, 1000)])")
30 call disp%show(
"if (0 /= getErrTableWrite(file = 'isFailedGeomCyclicFit.rnd.txt', table = rand)) error stop 'table output failed.'")
31 if (
0 /= getErrTableWrite(file
= 'isFailedGeomCyclicFit.rnd.txt', table
= rand))
error stop 'table output failed.'
32 call disp%show(
"call setUnique(getGeomCyclicRand(probSuccess, [(period, i = 1, 1000)]), unique = stepSuccess, count = freqSuccess, order = -1_IK)")
33 call setUnique(
getGeomCyclicRand(probSuccess, [(period, i
= 1,
1000)]), unique
= stepSuccess, count
= freqSuccess, order
= -1_IK)
35 call disp%show( stepSuccess ,
format = SK_
"(sp,20(g0,:,', '))")
37 call disp%show( freqSuccess ,
format = SK_
"(sp,20(g0,:,', '))")
38 call disp%show(
"if (isFailedGeomCyclicFit(stepSuccess, freqSuccess, period, probSuccessFit, normFacFit)) error stop 'Fitting failed.'")
39 if (
isFailedGeomCyclicFit(stepSuccess, freqSuccess, period, probSuccessFit, normFacFit))
error stop 'Fitting failed.'
40 call disp%show(
"[probSuccessFit, normFacFit]")
41 call disp%show( [probSuccessFit, normFacFit] )
42 call disp%show(
"x = getLinSpace(1._RKG, real(period, RKG), 500_IK) ! for visualization.")
44 call disp%show(
"y = normFacFit * probSuccessFit * (1 - probSuccessFit)**(x - 1) / (1 - (1 - probSuccessFit)**period) ! for visualization.")
45 y
= normFacFit
* probSuccessFit
* (
1 - probSuccessFit)
**(x
- 1)
/ (
1 - (
1 - probSuccessFit)
**period)
46 call disp%show(
"y = normFacFit * exp(getGeomCyclicLogPMF(x, probSuccessFit, period)) ! for visualization.")
48 call disp%show(
"if (0 /= getErrTableWrite(file = 'isFailedGeomCyclicFit.fit.txt', table = reshape([x, y], [size(x), 2]))) error stop 'table output failed.'")
49 if (
0 /= getErrTableWrite(file
= 'isFailedGeomCyclicFit.fit.txt', table
= reshape([x, y], [
size(x),
2])))
error stop 'table output failed.'
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Return a vector of unique values in the input array in place of the array itself.
Generate and return the natural logarithm of the Probability Mass Function (PMF) of the Cyclic Geomet...
Generate and return a scalar (or array of arbitrary rank of) random value(s) from the Cyclic Geometri...
Generate and return the iostat code resulting from writing the input table of rank 1 or 2 to the spec...
This is a generic method of the derived type display_type with pass attribute.
This is a generic method of the derived type display_type with pass attribute.
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
This module contains procedures and generic interfaces for finding unique values of an input array of...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Generate and return an object of type display_type.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
2probSuccess
= .
1; period
= 20
4if (
0 /= getErrTableWrite(file
= 'isFailedGeomCyclicFit.rnd.txt', table
= rand))
error stop 'table output failed.'
5call setUnique(
getGeomCyclicRand(probSuccess, [(period, i
= 1,
1000)]), unique
= stepSuccess, count
= freqSuccess, order
= -1_IK)
7+1,
+2,
+4,
+3,
+5,
+6,
+8,
+9,
+7,
+11,
+12,
+10,
+13,
+14,
+15,
+17,
+16,
+20,
+18,
+19
9+114,
+97,
+91,
+90,
+84,
+76,
+55,
+55,
+52,
+51,
+36,
+33,
+31,
+26,
+25,
+24,
+23,
+15,
+12,
+10
10if (
isFailedGeomCyclicFit(stepSuccess, freqSuccess, period, probSuccessFit, normFacFit))
error stop 'Fitting failed.'
11[probSuccessFit, normFacFit]
12+0.105651955925097061021975721637957888,
+976.572977671411732916580527403386634
14y
= normFacFit
* probSuccessFit
* (
1 - probSuccessFit)
**(x
- 1)
/ (
1 - (
1 - probSuccessFit)
**period)
16if (
0 /= getErrTableWrite(file
= 'isFailedGeomCyclicFit.fit.txt', table
= reshape([x, y], [
size(x),
2])))
error stop 'table output failed.'
Postprocessing of the example output ⛓
3import matplotlib.pyplot
as plt
16label = [
r"probSuccess = .2, period = 20"
20fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
23dfrnd = pd.read_csv(
"isFailedGeomCyclicFit.rnd.txt", delimiter =
",", header =
None)
24for j
in range(len(dfrnd.values[0,:])):
25 plt.hist( dfrnd.values[:,j]
26 , histtype =
"stepfilled"
31dffit = pd.read_csv(
"isFailedGeomCyclicFit.fit.txt", delimiter =
",", header =
None)
32for j
in range(1, len(dffit.values[0,:])):
33 plt.plot( dffit.values[:, 0]
41plt.xticks(fontsize = fontsize - 2)
42plt.yticks(fontsize = fontsize - 2)
44ax.set_xlabel(
"Cyclic Geometric Random Value ( integer-valued )", fontsize = 17)
45ax.set_ylabel(
"Count", fontsize = 17)
46ax.set_title(
"{} Cyclic Geometric random values and histogram fit".format(len(dfrnd.values[:, 0])), fontsize = 17)
48plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
49ax.tick_params(axis =
"y", which =
"minor")
50ax.tick_params(axis =
"x", which =
"minor")
52plt.savefig(
"isFailedGeomCyclicFit.png")
Visualization of the example output ⛓
- Test:
- test_pm_distGeomCyclic
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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- Author:
- Amir Shahmoradi, Monday March 6, 2017, 3:22 pm, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin.
Definition at line 2555 of file pm_distGeomCyclic.F90.