ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_distPois Module Reference

This module contains classes and procedures for computing various statistical quantities related to the Poisson distribution. More...

Data Types

type  distPois_type
 This is the derived type for signifying distributions that are of type Poisson as defined in the description of pm_distPois. More...
 
interface  getPoisCDF
 Generate and return the Cumulative Distribution Function (CDF) of the Poisson distribution for an input count within the discrete integer support of the distribution \([0, +\infty)\). More...
 
interface  getPoisLogPMF
 Generate and return the natural logarithm of the Probability Mass Function (PMF) of the Poisson distribution for an input count within the discrete integer support of the distribution \([0, +\infty)\). More...
 
interface  getPoisRand
 Generate and return a scalar (or array of arbitrary rank of) random value(s) from the Poisson distribution.
More...
 
interface  setPoisCDF
 Return the Cumulative Distribution Function (CDF) of the Poisson distribution. More...
 
interface  setPoisLogPMF
 Return the natural logarithm of the Probability Mass Function (PMF) of the Poisson distribution for an input count within the discrete integer support of the distribution \([0, +\infty)\). More...
 
interface  setPoisRand
 Return a scalar (or array of arbitrary rank of) random value(s) from the Poisson distribution.
More...
 

Variables

character(*, SK), parameter MODULE_NAME = "@pm_distPois"
 
real(RKB), parameter LAMBDA_LIMIT = 10._RKB
 The constant scalar of type real of kind RKB, representing the value of the parameter of the Poisson distribution above which the rejection method of Hormann, 1993, The transformed rejection method for generating Poisson random variables for generating Poisson-distributed random values is valid.
More...
 

Detailed Description

This module contains classes and procedures for computing various statistical quantities related to the Poisson distribution.

Specifically, this module contains routines for computing the following quantities of the Poisson distribution:

  1. the Probability Mass Function (PMF)
  2. the Cumulative Distribution Function (CDF)
  3. the Random Number Generation from the distribution (RNG)
  4. the Inverse Cumulative Distribution Function (ICDF) or the Quantile Function

The Poisson distribution is a discrete probability distribution that expresses the probability of a given number of events occurring in a fixed interval of time or space if these events occur with a known constant mean rate and independently of the time since the last event.
It is named after French mathematician Siméon Denis Poisson.
The Poisson distribution can also be used for the number of events in other specified interval types such as distance, area, or volume.
A discrete random variable \(X\) is said to have a Poisson distribution, with parameter \(\lambda > 0\) if it has a probability mass function given by,

\begin{equation} f(k; \lambda) = \pi(X = k) = \frac {\lambda^k \exp\left(-\lambda\right)}{k!} ~, \end{equation}

where

  1. \(k\) is the number of occurrences ( \(k = 0, 1, 2, \ldots\)), and
  2. \(\ms{!}\) is the factorial function.

The positive real number \(\lambda\) is equal to the expected value of \(X\) and also to its variance.

\begin{equation} \lambda = \up{E}(X) = \up{Var}(X) ~. \end{equation}

The CDF of the Poisson distribution with parameter \(\lambda\) is defined as,

\begin{eqnarray} \ms{CDF}(k | \lambda) &=& \frac{\Gamma(\lfloor k + 1 \rfloor, \lambda)}{\lfloor k \rfloor !} ~, \nonumber \\ &=& \exp\left(-\lambda\right) \sum _{j=0}^{\lfloor k \rfloor}{\frac{\lambda^{j}}{j!}} ~, \end{eqnarray}

where

  1. \(k\) is the number of occurrences ( \(k = 0, 1, 2, \ldots\)), and
  2. \(\ms{!}\) is the factorial function, and
  3. \(\Gamma(x, y) / \lfloor k \rfloor !\) is the regularized upper incomplete gamma function, and
  4. \(\lfloor k \rfloor\) is the floor function.

Random Number Generation

The RNG generic interfaces of this module use two different approaches for Poisson RNG for different ranges of \(\lambda\) parameter values of the Poisson distribution.

  1. When \(\lambda <\) LAMBDA_LIMIT, a RNG algorithm due to Donald Ervin Knuth is used.
  2. When LAMBDA_LIMIT \(\leq \lambda\), a rejection-based RNG algorithm due to Hormann, 1993, The transformed rejection method for generating Poisson random variables is used.
    This rejection method has an efficiency slight less than \(90\%\).
See also
pm_distPois
pm_distPois
Poisson distribution
Benchmarks:


Benchmark :: The runtime performance of getPoisLogPMF vs. setPoisLogPMF

1! Test the performance of `getPoisLogPMF()` vs. `setPoisLogPMF()`.
2program benchmark
3
4 use iso_fortran_env, only: error_unit
5 use pm_kind, only: SK, IK, RK, RKG => RK
7 use pm_bench, only: bench_type
8
9 implicit none
10
11 integer(IK) :: i
12 integer(IK) :: isize
13 integer(IK) :: fileUnit
14 integer(IK) , parameter :: NSIZE = 18_IK
15 integer(IK) , parameter :: NBENCH = 2_IK
16 integer(IK) :: arraySize(NSIZE)
17 integer(IK) , allocatable :: count(:)
18 real(RKG) , allocatable :: array(:)
19 real(RKG) :: dummy = 0._RKG
20 type(bench_type) :: bench(NBENCH)
21 type(xoshiro256ssw_type) :: rng
22
23 rng = xoshiro256ssw_type()
24
25 bench(1) = bench_type(name = SK_"getPoisLogPMF", exec = getPoisLogPMF , overhead = setOverhead)
26 bench(2) = bench_type(name = SK_"setPoisLogPMF", exec = setPoisLogPMF , overhead = setOverhead)
27
28 arraySize = [( 2_IK**isize, isize = 1_IK, NSIZE )]
29
30 write(*,"(*(g0,:,' '))")
31 write(*,"(*(g0,:,' vs. '))") (bench(i)%name, i = 1, NBENCH)
32 write(*,"(*(g0,:,' '))")
33
34 open(newunit = fileUnit, file = "main.out", status = "replace")
35
36 write(fileUnit, "(*(g0,:,','))") "arraySize", (bench(i)%name, i = 1, NBENCH)
37
38 loopOverarraySize: do isize = 1, NSIZE
39
40 write(*,"(*(g0,:,' '))") "Benchmarking with size", arraySize(isize)
41
42 allocate(array(arraySize(isize)), count(arraySize(isize)))
43 do i = 1, NBENCH
44 bench(i)%timing = bench(i)%getTiming(minsec = 0.05_RK)
45 end do
46 deallocate(array, count)
47
48 write(fileUnit,"(*(g0,:,','))") arraySize(isize), (bench(i)%timing%mean, i = 1, NBENCH)
49
50 end do loopOverarraySize
51
52 write(*,"(*(g0,:,' '))") dummy
53 write(*,"(*(g0,:,' '))")
54
55 close(fileUnit)
56
57contains
58
59 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
60 ! procedure wrappers.
61 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
62
63 subroutine setOverhead()
64 call initialize()
65 call finalize()
66 end subroutine
67
68 subroutine initialize()
69 use pm_distUnif, only: setUnifRand
70 call setUnifRand(rng, count, 0_IK, 1023_IK)
71 end subroutine
72
73 subroutine finalize()
74 dummy = dummy + array(1)
75 end subroutine
76
77 subroutine setPoisLogPMF()
78 block
79 use pm_distPois, only: setPoisLogPMF
80 call initialize()
81 call setPoisLogPMF(array, count, lambda = 1._RKG)
82 call finalize()
83 end block
84 end subroutine
85
86 subroutine getPoisLogPMF()
87 block
88 use pm_distPois, only: getPoisLogPMF
89 call initialize()
90 array = getPoisLogPMF(count, lambda = 1._RKG)
91 call finalize()
92 end block
93 end subroutine
94
95end program benchmark
Generate and return an object of type timing_type containing the benchmark timing information and sta...
Definition: pm_bench.F90:574
Generate and return the natural logarithm of the Probability Mass Function (PMF) of the Poisson distr...
Return the natural logarithm of the Probability Mass Function (PMF) of the Poisson distribution for a...
Return a uniform random scalar or contiguous array of arbitrary rank of randomly uniformly distribute...
This module contains abstract interfaces and types that facilitate benchmarking of different procedur...
Definition: pm_bench.F90:41
This module contains classes and procedures for computing various statistical quantities related to t...
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...
Definition: pm_kind.F90:268
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
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
This is the class for creating benchmark and performance-profiling objects.
Definition: pm_bench.F90:386
This is the derived type for declaring and generating objects of type xoshiro256ssw_type containing a...

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

Postprocessing of the benchmark output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6
7fontsize = 14
8
9methods = ["setPoisLogPMF", "getPoisLogPMF"]
10
11df = pd.read_csv("main.out")
12
13
16
17ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
18ax = plt.subplot()
19
20for method in methods:
21 plt.plot( df["arraySize"].values
22 , df[method].values
23 , linewidth = 2
24 )
25
26plt.xticks(fontsize = fontsize)
27plt.yticks(fontsize = fontsize)
28ax.set_xlabel("Array Size", fontsize = fontsize)
29ax.set_ylabel("Runtime [ seconds ]", fontsize = fontsize)
30ax.set_title("setPoisLogPMF() vs. getPoisLogPMF()\nLower is better.", fontsize = fontsize)
31ax.set_xscale("log")
32ax.set_yscale("log")
33plt.minorticks_on()
34plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
35ax.tick_params(axis = "y", which = "minor")
36ax.tick_params(axis = "x", which = "minor")
37ax.legend ( methods
38 #, loc='center left'
39 #, bbox_to_anchor=(1, 0.5)
40 , fontsize = fontsize
41 )
42
43plt.tight_layout()
44plt.savefig("benchmark.getPoisLogPMF_vs_setPoisLogPMF.runtime.png")
45
46
49
50ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
51ax = plt.subplot()
52
53plt.plot( df["arraySize"].values
54 , np.ones(len(df["arraySize"].values))
55 #, linestyle = "--"
56 #, color = "black"
57 , linewidth = 2
58 )
59plt.plot( df["arraySize"].values
60 , df["getPoisLogPMF"].values / df["setPoisLogPMF"].values
61 , linewidth = 2
62 )
63
64plt.xticks(fontsize = fontsize)
65plt.yticks(fontsize = fontsize)
66ax.set_xlabel("Array Size", fontsize = fontsize)
67ax.set_ylabel("Runtime compared to setPoisLogPMF()", fontsize = fontsize)
68ax.set_title("getPoisLogPMF() / setPoisLogPMF()\nLower means faster. Lower than 1 means faster than setPoisLogPMF().", fontsize = fontsize)
69ax.set_xscale("log")
70#ax.set_yscale("log")
71plt.minorticks_on()
72plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
73ax.tick_params(axis = "y", which = "minor")
74ax.tick_params(axis = "x", which = "minor")
75ax.legend ( ["setPoisLogPMF", "getPoisLogPMF"]
76 #, bbox_to_anchor = (1, 0.5)
77 #, loc = "center left"
78 , fontsize = fontsize
79 )
80
81plt.tight_layout()
82plt.savefig("benchmark.getPoisLogPMF_vs_setPoisLogPMF.runtime.ratio.png")

Visualization of the benchmark output

Benchmark moral
  1. The procedures under the generic interface getPoisLogPMF are functions while the procedures under the generic interface setPoisLogPMF are subroutines.
    From the benchmark results, it appears that the functional interface performs slightly less efficiently than the subroutine interface when the input array size is small.
    Otherwise, the difference appears to be marginal and insignificant in most practical situations.


Benchmark :: The runtime performance of setPoisLogPMF with and without logLambda

1program benchmark
2
3 use iso_fortran_env, only: error_unit
4 use pm_bench, only: bench_type
6 use pm_distUnif, only: setUnifRand
7 use pm_kind, only: SK, IK, RK, RKG => RK
8
9 implicit none
10
11 integer(IK) :: i
12 integer(IK) :: isize
13 integer(IK) :: fileUnit
14 integer(IK) , parameter :: NSIZE = 18_IK
15 integer(IK) , parameter :: NBENCH = 2_IK
16 integer(IK) :: arraySize(0:NSIZE)
17 integer(IK) , allocatable :: count(:)
18 real(RKG) , allocatable :: logPMF(:)
19 real(RKG) , allocatable :: logLambda(:), lambda(:)
20 real(RKG) :: dummy = 0._RKG
21 type(bench_type) :: bench(NBENCH)
22 type(xoshiro256ssw_type) :: rng
23
24 rng = xoshiro256ssw_type()
25
26 bench(1) = bench_type(name = SK_"logLambdaMissing", exec = logLambdaMissing , overhead = setOverhead)
27 bench(2) = bench_type(name = SK_"logLambdaPresent", exec = logLambdaPresent , overhead = setOverhead)
28
29 arraySize = [( 2_IK**isize, isize = 0_IK, NSIZE )]
30
31 write(*,"(*(g0,:,' '))")
32 write(*,"(*(g0,:,' vs. '))") (bench(i)%name, i = 1, NBENCH)
33 write(*,"(*(g0,:,' '))")
34
35 open(newunit = fileUnit, file = "main.out", status = "replace")
36
37 write(fileUnit, "(*(g0,:,','))") "arraySize", (bench(i)%name, i = 1, NBENCH)
38
39 loopOverarraySize: do isize = 0, NSIZE
40
41 write(*,"(*(g0,:,' '))") "Benchmarking with size", arraySize(isize)
42
43 allocate(logPMF(arraySize(isize)), count(arraySize(isize)), lambda(arraySize(isize)))
44 call setUnifRand(rng, count, 0_IK, 1023_IK)
45 call random_number(lambda)
46 lambda = 1._RKG - lambda
47 logLambda = log(lambda)
48 do i = 1, NBENCH
49 bench(i)%timing = bench(i)%getTiming(minsec = 0.05_RK)
50 end do
51 deallocate(logPMF, count, lambda)
52
53 write(fileUnit,"(*(g0,:,','))") arraySize(isize), (bench(i)%timing%mean, i = 1, NBENCH)
54
55 end do loopOverarraySize
56
57 write(*,"(*(g0,:,' '))") dummy
58 write(*,"(*(g0,:,' '))")
59
60 close(fileUnit)
61
62contains
63
64 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
65 ! procedure wrappers.
66 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
67
68 subroutine setOverhead()
69 call initialize()
70 call finalize()
71 end subroutine
72
73 subroutine initialize()
74 !call random_number(lambda)
75 !logLambda = log(lambda)
76 end subroutine
77
78 subroutine finalize()
79 dummy = dummy + logPMF(1)
80 end subroutine
81
82 subroutine logLambdaPresent()
83 use pm_distPois, only: setPoisLogPMF
84 call initialize()
85 if (arraySize(isize) > 1_IK) then
86 call setPoisLogPMF(logPMF, count, lambda, logLambda)
87 else
88 call setPoisLogPMF(logPMF(1), count(1), lambda(1), logLambda(1))
89 end if
90 call finalize()
91 end subroutine
92
93 subroutine logLambdaMissing()
94 use pm_distPois, only: getPoisLogPMF
95 call initialize()
96 if (arraySize(isize) > 1_IK) then
97 logPMF = getPoisLogPMF(count, lambda = lambda)
98 else
99 logPMF(1) = getPoisLogPMF(count(1), lambda = lambda(1))
100 end if
101 call finalize()
102 end subroutine
103
104end program benchmark

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

Postprocessing of the benchmark output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6
7fontsize = 14
8
9methods = ["logLambdaPresent", "logLambdaMissing"]
10
11df = pd.read_csv("main.out")
12
13
16
17ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
18ax = plt.subplot()
19
20for method in methods:
21 plt.plot( df["arraySize"].values
22 , df[method].values
23 , linewidth = 2
24 )
25
26plt.xticks(fontsize = fontsize)
27plt.yticks(fontsize = fontsize)
28ax.set_xlabel("Array Size", fontsize = fontsize)
29ax.set_ylabel("Runtime [ seconds ]", fontsize = fontsize)
30ax.set_title("logLambdaPresent() vs. logLambdaMissing()\nLower is better.", fontsize = fontsize)
31ax.set_xscale("log")
32ax.set_yscale("log")
33plt.minorticks_on()
34plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
35ax.tick_params(axis = "y", which = "minor")
36ax.tick_params(axis = "x", which = "minor")
37ax.legend ( methods
38 #, loc='center left'
39 #, bbox_to_anchor=(1, 0.5)
40 , fontsize = fontsize
41 )
42
43plt.tight_layout()
44plt.savefig("benchmark.setPoisLogPMF-logLambda-missing_vs_present.runtime.png")
45
46
49
50ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
51ax = plt.subplot()
52
53plt.plot( df["arraySize"].values
54 , np.ones(len(df["arraySize"].values))
55 #, linestyle = "--"
56 #, color = "black"
57 , linewidth = 2
58 )
59plt.plot( df["arraySize"].values
60 , df["logLambdaMissing"].values / df["logLambdaPresent"].values
61 , linewidth = 2
62 )
63
64plt.xticks(fontsize = fontsize)
65plt.yticks(fontsize = fontsize)
66ax.set_xlabel("Array Size", fontsize = fontsize)
67ax.set_ylabel("Runtime compared to logLambdaPresent()", fontsize = fontsize)
68ax.set_title("logLambdaMissing / logLambdaPresent\nLower means faster. Lower than 1 means faster than logLambdaPresent.", fontsize = fontsize)
69ax.set_xscale("log")
70#ax.set_yscale("log")
71plt.minorticks_on()
72plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
73ax.tick_params(axis = "y", which = "minor")
74ax.tick_params(axis = "x", which = "minor")
75ax.legend ( ["logLambdaPresent", "logLambdaMissing"]
76 #, bbox_to_anchor = (1, 0.5)
77 #, loc = "center left"
78 , fontsize = fontsize
79 )
80
81plt.tight_layout()
82plt.savefig("benchmark.setPoisLogPMF-logLambda-missing_vs_present.runtime.ratio.png")

Visualization of the benchmark output

Benchmark moral
  1. The procedures under the generic interface setPoisLogPMF accept an extra argument logLambda = log(lambda) while the procedures under the generic interface getPoisLogPMF compute this term internally with every procedure call.
    In the presence of this argument, the logarithmic computation log(lambda) will be avoided.
    As such, the presence of logLambda is expected to lead to faster computations.
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.

  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

Variable Documentation

◆ LAMBDA_LIMIT

real(RKB), parameter pm_distPois::LAMBDA_LIMIT = 10._RKB

The constant scalar of type real of kind RKB, representing the value of the parameter of the Poisson distribution above which the rejection method of Hormann, 1993, The transformed rejection method for generating Poisson random variables for generating Poisson-distributed random values is valid.

This constant exists merely as a reference with which decision can be made about the proper setPoisRand interface usage.

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.

  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 762 of file pm_distPois.F90.

◆ MODULE_NAME

character(*, SK), parameter pm_distPois::MODULE_NAME = "@pm_distPois"

Definition at line 129 of file pm_distPois.F90.