Generate and return a scalar (or array of arbitrary rank of) random value(s) from the Poisson distribution.
More...
Generate and return a scalar (or array of arbitrary rank of) random value(s) from the Poisson distribution.
See the documentation of pm_distPois for more details.
- Parameters
-
[in] | lambda | : The input scalar (or array of the same rank, shape, and size as other array-like arguments), of
-
type
real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
representing the location/scale parameter of the distribution.
|
- Returns
rand
: The output positive scalar (or array of the same rank, shape, and size as other array-like arguments), of type integer
of default kind IK, containing random value(s) from the distribution.
Possible calling interfaces ⛓
integer(IK) :: rand
Generate and return a scalar (or array of arbitrary rank of) random value(s) from the Poisson distrib...
This module contains classes and procedures for computing various statistical quantities related to t...
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...
- Warning
- The condition
0 < lambda
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
- setPoisRand
Example usage ⛓
12 integer(IK),
parameter :: NP
= 1000_IK
13 integer(IK) :: rand(NP)
14 real(RKG) :: lambda(NP)
16 type(display_type) :: disp
19 call setLinSpace(lambda, x1
= 0.1_RKG, x2
= 100._RKG)
24 call disp%show(
"rand(1) = getPoisRand(lambda(1))")
33 call disp%show(
"rand(1:2) = getPoisRand(lambda(1))")
42 call disp%show(
"rand(1:2) = getPoisRand(lambda(1))")
51 call disp%show(
"rand(NP-2:NP) = getPoisRand(lambda(NP-2:NP))")
62 integer(IK) :: fileUnit, i
63 integer(IK) ,
parameter :: NP
= 5000_IK
64 real(RKG) ,
parameter :: lambda(
4)
= [.
1_RKG,
1._RKG,
4._RKG,
11._RKG]
65 open(newunit
= fileUnit, file
= "getPoisRand.IK.txt")
Return the linSpace output argument with size(linSpace) elements of evenly-spaced values over the int...
Return the logSpace output argument with size(logSpace) elements of logarithmically-evenly-spaced val...
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 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.
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
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 ⛓
24+99.8000031,
+99.9000015,
+100.000000
Postprocessing of the example output ⛓
3import matplotlib.pyplot
as plt
16xlab = {
"CK" :
"Poisson Random Value ( real/imaginary components )"
17 ,
"IK" :
"Poisson Random Value ( integer-valued )"
18 ,
"RK" :
"Poisson Random Value ( real-valued )"
20legends = [
r"$\lambda = 0.1$"
26for kind
in [
"IK",
"CK",
"RK"]:
28 pattern =
"*." + kind +
".txt"
29 fileList = glob.glob(pattern)
30 if len(fileList) == 1:
32 df = pd.read_csv(fileList[0], delimiter =
" ", header =
None)
34 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
37 for j
in range(len(df.values[0,:])):
39 plt.hist( df.values[:,j]
40 , histtype =
"stepfilled"
45 plt.hist( df.values[:,j]
46 , histtype =
"stepfilled"
53 plt.xticks(fontsize = fontsize - 2)
54 plt.yticks(fontsize = fontsize - 2)
55 ax.set_xlabel(xlab[kind], fontsize = 17)
56 ax.set_ylabel(
"Count", fontsize = 17)
57 ax.set_title(
"Histograms of {} Poisson random values".format(len(df.values[:, 0])), fontsize = 17)
59 plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
60 ax.tick_params(axis =
"y", which =
"minor")
61 ax.tick_params(axis =
"x", which =
"minor")
63 plt.savefig(fileList[0].replace(
".txt",
".png"))
65 elif len(fileList) > 1:
67 sys.exit(
"Ambiguous file list exists.")
Visualization of the example output ⛓
- Test:
- test_pm_distPois
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, Oct 16, 2009, 11:14 AM, Michigan
Definition at line 821 of file pm_distPois.F90.