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

Generate and return the natural logarithm of an approximation of the Probability Density Function (PDF) of the Multiple MultiVariate Uniformly Ellipsoidal (MMVUE) Distribution.
More...

Detailed Description

Generate and return the natural logarithm of an approximation of the Probability Density Function (PDF) of the Multiple MultiVariate Uniformly Ellipsoidal (MMVUE) Distribution.

See the documentation of pm_distUnifElls for details of the definition of the PDF.

Note
Note that the procedures of this generic interface do not require an input X value representing a location within the domain of the density function.
This is fine and intentional by design because the PDF is uniform across the entire support of the PDF.
Parameters
[in,out]rng: The input/output scalar that can be an object of,
  1. type rngf_type, implying the use of intrinsic Fortran uniform RNG.
  2. type xoshiro256ssw_type, implying the use of xoshiro256** uniform RNG.
[in]mean: The input contiguous matrix of shape (1:ndim, 1:nell), of the same type and kind as the output logPDF, representing the centers of the ellipsoids in the distribution.
[in]chol: The input contiguous array of shape (1:ndim, 1:ndim, 1:nell) whose specified triangular subset contains the Cholesky Factorization of the Gramian matrix of the MMVUE distribution.
[in]subset: The input scalar constant that can be any of the following:
  1. the constant uppDia or an object of type uppDia_type implying that the upper-diagonal triangular block of the input chol must be used while the lower subset is not referenced.
  2. the constant lowDia or an object of type lowDia_type implying that the lower-diagonal triangular block of the input chol must be used while the upper subset is not referenced.
This argument is merely a convenience to differentiate the different procedure functionalities within this generic interface.
[in]invGram: The input array of shape (1:ndim, 1:ndim, 1:nell) of the same type and kind as the input argument point, containing the collection of square matrices of full inverse representative Gramian matrix of the \(\ndim\)-dimensional ellipsoids in the distribution.
This argument is needed to determine the membership of the generated random vectors for estimating the effective volumes the ellipsoids.
[in]nsim: The input positive scalar of type integer of default kind IK, containing the number of random number simulations in approximating the PDF of the distribution.
The larger this integer, the more accurate the logPDF estimate will be.
(optional, default = 10000.)
[in]normed: The input positive scalar of type logical of default kind LK.
  1. If set to .false., the output logPDF will contain the usual normalizing factor corresponding to the volume of a unit-ball in the corresponding number of dimensions of the domain of the distribution.
  2. If set to .true., the output logPDF will not contain the factor corresponding to the volume of a unit-ball in the corresponding number of dimensions of the domain of the distribution.
    This factor is returned by the generic interface getLogVolUnitBall.
    In such a case, the output logPDF is not simply the natural logarithm of the inverse of effective volumes of the ellipsoids, but instead, their effective sum of the square-roots of the multiplicative traces of their Cholesky factors.
    This normalized value is occasionally useful in clustering algorithms or when the normalization factor is irrelevant.
    Beware that in this case, the output logPDF does not represent the actual log(PDF) of the distribution, but is a constant-added equivalent to it.
(optional, default = .false..)
Returns
logPDF : The output scalar of,
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the natural logarithm of the PDF of the distribution.


Possible calling interfaces

logPDF = getUnifEllsLogPDF(rng, mean(1:ndim, 1:nell), nsim = nsim, normed = normed) ! all ellipsoids are unit-balls centered on `mean` slices.
logPDF = getUnifEllsLogPDF(rng, mean(1:ndim, 1:nell), chol(1:ndim, 1:ndim, 1:nell), subset, invGram(1:ndim, 1:ndim, 1:nell), nsim = nsim, normed = normed)
!
Generate and return the natural logarithm of an approximation of the Probability Density Function (PD...
This module contains classes and procedures for computing various statistical quantities related to t...
Warning
The condition 0 < nsim must hold for the corresponding input arguments.
The condition all(shape(invGram) == [size(mean, 1), size(mean, 1), size(mean, 2)) must hold for the corresponding input arguments.
The condition size(chol, 1) <= size(chol, 2) must hold for the corresponding input arguments (to ensure that arguments can be passed without data copy).
The condition all([size(chol, 1), size(chol, 3)] == shape(mean)) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The procedures under discussion are impure. The procedures of this generic interface are pure when the input argument rng is set to xoshiro256ssw_type and the compile-time macro CHECK_ENABLED is set to 0 or is undefined.
Warning
The condition 0 < nell must hold for the corresponding input arguments.
The condition all(invmul <= 1) must hold for the corresponding input arguments.
The condition size(chol, 1) == size(chol, 2) 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
pm_ellipsoid
setUnifEllsRand


Example usage

1program example
2
3 use pm_kind, only: SK
4 use pm_kind, only: IK, LK
5 use pm_io, only: display_type
8 use pm_matrixChol, only: getMatChol, uppDia, lowDia
10 use pm_distChol, only: getCholRand, uppDia
11 use pm_matrixInv, only: getMatInv, choUpp
12 use pm_io, only: getErrTableWrite, trans
13
14 implicit none
15
16 integer(IK) :: ndim, iell, nell, nsim
17
18 type(display_type) :: disp
19 disp = display_type(file = "main.out.F90")
20
21 block
22 use pm_kind, only: TKG => RKD
23 real(TKG), allocatable :: logPDF, mean(:,:), chol(:, :, :), invGram(:, :, :)
24
25 call disp%skip()
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show("! Random vectors from Multiple Multivariate Uniform Balls with particular means.")
28 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 call disp%skip()
30
31 call disp%skip()
32 call disp%show("ndim = 2; nell = getUnifRand(5, 10); nsim = 2000")
33 ndim = 2; nell = getUnifRand(5, 10); nsim = 2000
34 call disp%show("[ndim, nell, nsim]")
35 call disp%show( [ndim, nell, nsim] )
36 call disp%show("mean = getUnifRand(-5._TKG, +5._TKG, ndim, nell)")
37 mean = getUnifRand(-5._TKG, +5._TKG, ndim, nell)
38 call disp%show("mean")
39 call disp%show( mean )
40 call disp%show("logPDF = getUnifEllsLogPDF(rngf_type(), mean)")
41 logPDF = getUnifEllsLogPDF(rngf_type(), mean)
42 call disp%show("logPDF")
43 call disp%show( logPDF )
44 call disp%skip()
45
46 call disp%skip()
47 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
48 call disp%show("! Random vectors from Multiple Multivariate Uniform Ellipsoids with particular means.")
49 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
50 call disp%skip()
51
52 call disp%show("call setResized(chol, [ndim, ndim, nell])")
53 call setResized(chol, [ndim, ndim, nell])
54 call disp%show("call setResized(invGram, [ndim, ndim, nell])")
55 call setResized(invGram, [ndim, ndim, nell])
56 call disp%show("do iell = 1, nell; chol(:, :, iell) = getCholRand(0._TKG, ndim, uppDia); invGram(:, :, iell) = getMatInv(chol(:, :, iell), choUpp); end do")
57 do iell = 1, nell; chol(:, :, iell) = getCholRand(0._TKG, ndim, uppDia); invGram(:, :, iell) = getMatInv(chol(:, :, iell), choUpp); end do
58 call disp%show("logPDF = getUnifEllsLogPDF(rngf_type(), mean, chol, uppDia, invGram)")
59 logPDF = getUnifEllsLogPDF(rngf_type(), mean, chol, uppDia, invGram)
60 call disp%show("logPDF")
61 call disp%show( logPDF )
62 call disp%skip()
63
64 end block
65
66end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate and return a random upper and lower Cholesky factorization.
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
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
Generate and return the upper or the lower Cholesky factorization of the input symmetric positive-def...
Generate and return the full inverse of an input matrix of general or triangular form directly or thr...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains classes and procedures for generating random upper or lower Cholesky factor tria...
Definition: pm_distChol.F90:43
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 RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
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 module contains procedures and generic interfaces for computing the Cholesky factorization of po...
This module contains abstract and concrete derived types and procedures related to the inversion of s...
This is a concrete derived type whose instances can be used to define/request the default uniform ran...
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! Random vectors from Multiple Multivariate Uniform Balls with particular means.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7ndim = 2; nell = getUnifRand(5, 10); nsim = 2000
8[ndim, nell, nsim]
9+2, +6, +2000
10mean = getUnifRand(-5._TKG, +5._TKG, ndim, nell)
11mean
12+2.9228221659175375, +0.46951155098910924, -2.1250736221240052, +2.3044603380669448, -3.6604889981740367, +3.4387854709374786
13-0.85015055710606369, +0.83266500477845096, -4.8240433032172731, +0.97730118318714254, +2.4205490950935937, -0.17029210514103088
14logPDF = getUnifEllsLogPDF(rngf_type(), mean)
15logPDF
16+1.7326992480985099
17
18
19!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
20! Random vectors from Multiple Multivariate Uniform Ellipsoids with particular means.
21!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
22
23call setResized(chol, [ndim, ndim, nell])
24call setResized(invGram, [ndim, ndim, nell])
25do iell = 1, nell; chol(:, :, iell) = getCholRand(0._TKG, ndim, uppDia); invGram(:, :, iell) = getMatInv(chol(:, :, iell), choUpp); end do
26logPDF = getUnifEllsLogPDF(rngf_type(), mean, chol, uppDia, invGram)
27logPDF
28+2.4145904206339979
29
30
Test:
test_pm_distUnifElls


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, Monday March 6, 2017, 3:22 pm, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin.

Definition at line 310 of file pm_distUnifElls.F90.


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