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

Return the natural logarithm of the volume of an \(\ndim\)-dimensional ball of unit-radius. More...

Detailed Description

Return the natural logarithm of the volume of an \(\ndim\)-dimensional ball of unit-radius.

The computation of the volume of an n-ball requires involves factorials which are computed in the procedures of this generic interface via the Fortran intrinsic log_gamma().
This generic subroutine interface is an exact functional-interface replacement for the generic functional interface getLogVolUnitBall.

Parameters
[out]logVolUnitBall: The output scalar (or array of the same rank as other input array-like arguments) of,
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing natural logarithm of the volume of the unit-radius hyper-ball.
[in]ndim: The input scalar containing the number of dimensions of the unit-radius hyper-ball.
It can be,
  1. of type integer of default kind IK, in which case the output will be computed using an iterative factorial algorithm, or
  2. of the same type and kind as the output logVolUnitBall, in which case the output will be computed using the Fortran intrinsic log_gamma().


Possible calling interfaces

call setLogVolUnitBall(logVolUnitBall, ndim)
Return the natural logarithm of the volume of an -dimensional ball of unit-radius.
This module contains classes and procedures for setting up and computing the properties of the hyper-...
Warning
The condition 0 <= ndim must hold for the corresponding input arguments.
This condition is 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.
Remarks
The procedures under discussion are elemental.
See also
getLogVolUnitBall
Volume of an n-ball
Particular values of the Gamma function


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKG => RKH
4 use pm_io, only: display_type
6
7 implicit none
8
9 real(RKG) :: logVolUnitBall(5)
10
11 type(display_type) :: disp
12 disp = display_type(file = "main.out.F90")
13
14 call disp%skip()
15 call disp%show("call setLogVolUnitBall(logVolUnitBall(1), ndim = 2_IK)")
16 call setLogVolUnitBall(logVolUnitBall(1), ndim = 2_IK)
17 call disp%show("logVolUnitBall(1)")
18 call disp%show( logVolUnitBall(1) )
19 call disp%skip()
20
21 call disp%skip()
22 call disp%show("call setLogVolUnitBall(logVolUnitBall(1), ndim = 2._RKG)")
23 call setLogVolUnitBall(logVolUnitBall(1), ndim = 2._RKG)
24 call disp%show("logVolUnitBall(1)")
25 call disp%show( logVolUnitBall(1) )
26 call disp%skip()
27
28 call disp%skip()
29 call disp%show("call setLogVolUnitBall(logVolUnitBall(1:5), ndim = [integer(IK) :: 0, 1, 2, 3, 4])")
30 call setLogVolUnitBall(logVolUnitBall(1:5), ndim = [integer(IK) :: 0, 1, 2, 3, 4])
31 call disp%show("logVolUnitBall(1:5)")
32 call disp%show( logVolUnitBall(1:5) )
33 call disp%skip()
34
35 call disp%skip()
36 call disp%show("call setLogVolUnitBall(logVolUnitBall(1:5), ndim = [real(RKG) :: 0, 1, 2, 3, 4.5])")
37 call setLogVolUnitBall(logVolUnitBall(1:5), ndim = [real(RKG) :: 0, 1, 2, 3, 4.5])
38 call disp%show("logVolUnitBall(1:5)")
39 call disp%show( logVolUnitBall(1:5) )
40 call disp%skip()
41
42 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
43 ! Output an example array for visualization.
44 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
45
46 block
47 use pm_arraySpace, only: setLinSpace
48 integer(IK), parameter :: nsim = 1000_IK
49 integer(IK) :: fileUnit, i
50 real :: volUnitBall(nsim), ndim(nsim)
51 open(newunit = fileUnit, file = "setLogVolUnitBall.RK.txt")
52 call setLinSpace(ndim, 0., 20.)
53 call setLogVolUnitBall(volUnitBall, ndim)
54 do i = 1, nsim
55 write(fileUnit, "(*(g0,:,','))") ndim(i), exp(volUnitBall(i))
56 end do
57 close(fileUnit)
58 end block
59
60end program example
Return the linSpace output argument with size(linSpace) elements of evenly-spaced values over the int...
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 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 RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
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
2call setLogVolUnitBall(logVolUnitBall(1), ndim = 2_IK)
3logVolUnitBall(1)
4+1.14472988584940017414342735135305866
5
6
7call setLogVolUnitBall(logVolUnitBall(1), ndim = 2._RKG)
8logVolUnitBall(1)
9+1.14472988584940017414342735135305866
10
11
12call setLogVolUnitBall(logVolUnitBall(1:5), ndim = [integer(IK) :: 0, 1, 2, 3, 4])
13logVolUnitBall(1:5)
14+0.00000000000000000000000000000000000, +0.693147180559945309417232121458176575, +1.14472988584940017414342735135305866, +1.43241195830118110158264635734688620, +1.59631259113885503886962258124794093
15
16
17call setLogVolUnitBall(logVolUnitBall(1:5), ndim = [real(RKG) :: 0, 1, 2, 3, 4.5])
18logVolUnitBall(1:5)
19+0.00000000000000000000000000000000000, +0.693147180559945309417232121458176575, +1.14472988584940017414342735135305866, +1.43241195830118110158264635734688600, +1.63984031205242503356424402200248520
20
21

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 ["RK"]:
13
14 pattern = "*." + kind + ".txt"
15 fileList = glob.glob(pattern)
16
17 for file in fileList:
18
19 df = pd.read_csv(file, delimiter = ",", header = None)
20
21 # definitions for the axes
22 #left, width = 0.1, 0.65
23 #bottom, height = 0.1, 0.65
24 #spacing = 0.015
25
26 # start with a square Figure
27 fig = plt.figure() # figsize = (8, 6)
28
29 plt.rcParams.update({'font.size': fontsize - 2})
30 ax = plt.subplot()
31 #ax = fig.add_axes([left, bottom, width, height]) # scatter plot
32 #ax_histx = fig.add_axes([left, bottom + height + spacing, width, 0.2], sharex = ax) # histx
33 #ax_histy = fig.add_axes([left + width + spacing, bottom, 0.2, height], sharey = ax) # histy
34 #
35 #for axes in [ax, ax_histx, ax_histy]:
36 # axes.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
37 # axes.tick_params(axis = "y", which = "minor")
38 # axes.tick_params(axis = "x", which = "minor")
39 #
40
43
44 # the scatter plot:
45 for i in range(0, len(df.values[0,:]), 2):
46 ax.scatter ( df.values[:,i]
47 , df.values[:,i+1]
48 , s = 8
49 , zorder = 1000
50 )
51
52 #ax_histx.hist(df.values[:, 0], bins = 50, zorder = 1000)
53 #ax_histy.hist(df.values[:, 1], bins = 50, orientation = "horizontal", zorder = 1000)
54
55 #ax.set_aspect("equal")
56 ax.set_xlim([0, df.values[-1,:].max()])
57 ax.set_xlabel("ndim", fontsize = 17)
58 ax.set_ylabel("Volume of ndim unit-ball", fontsize = 17)
59 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
60 plt.minorticks_on()
61
62 plt.tight_layout()
63 plt.savefig(file.replace(".txt",".png"))

Visualization of the example output
Benchmarks:


Benchmark :: The runtime performance of setLogVolUnitBall for integer vs. real input ndim.

1! Test the performance of `setLogVolUnitBall()` with and without the selection `control` argument.
2program benchmark
3
4 use pm_bench, only: bench_type
6 use pm_kind, only: IK, LK, RK1 => RKD, RK2 => RKS, RK, SK
7 use iso_fortran_env, only: error_unit
8
9 implicit none
10
11 integer(IK) :: i
12 integer(IK) :: idim
13 integer(IK) :: fileUnit
14 integer(IK) :: ndim_IK
15
16 real(RK1) :: logVolUnitBall_RK1
17 real(RK1) :: dum_RK1 = 0._RK1
18 real(RK1) :: ndim_RK1
19
20 real(RK2) :: logVolUnitBall_RK2
21 real(RK2) :: dum_RK2 = 0._RK2
22 real(RK2) :: ndim_RK2
23 type(bench_type) , allocatable :: bench(:)
24
25 bench = [ bench_type(name = SK_"ndimInt2RK1", exec = ndimInt2RK1, overhead = setOverhead_RK1) &
26 , bench_type(name = SK_"ndimInt2RK2", exec = ndimInt2RK2, overhead = setOverhead_RK1) &
27 , bench_type(name = SK_"ndimRealRK1", exec = ndimRealRK1, overhead = setOverhead_RK1) &
28 , bench_type(name = SK_"ndimRealRK2", exec = ndimRealRK2, overhead = setOverhead_RK2) &
29 ]
30
31 write(*,"(*(g0,:,' '))")
32 write(*,"(*(g0,:,' '))") "setLogVolUnitBallInt() vs. setLogVolUnitBallReal()"
33 write(*,"(*(g0,:,' '))")
34
35 open(newunit = fileUnit, file = "main.out", status = "replace")
36
37 write(fileUnit, "(*(g0,:,','))") "ndim", (bench(i)%name, i = 1, size(bench))
38 loopOverArraySize: do idim = 0, 20_IK
39
40 ndim_IK = idim
41 ndim_RK1 = real(idim, RK1)
42 ndim_RK2 = real(idim, RK2)
43 write(*,"(*(g0,:,' '))") "Benchmarking with ndim = ", idim
44 do i = 1, size(bench)
45 bench(i)%timing = bench(i)%getTiming(minsec = 0.04_RK)
46 end do
47 write(fileUnit,"(*(g0,:,','))") idim, (bench(i)%timing%mean, i = 1, size(bench))
48
49 end do loopOverArraySize
50 write(*,"(*(g0,:,' '))") dum_RK1
51 write(*,"(*(g0,:,' '))")
52
53 close(fileUnit)
54
55contains
56
57 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
58 ! procedure wrappers.
59 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
60
61 subroutine setOverhead_RK1()
62 call getDummy_RK1()
63 end subroutine
64
65 subroutine setOverhead_RK2()
66 call getDummy_RK2()
67 end subroutine
68
69 subroutine getDummy_RK1()
70 dum_RK1 = dum_RK1 + logVolUnitBall_RK1
71 end subroutine
72
73 subroutine getDummy_RK2()
74 dum_RK2 = dum_RK2 + logVolUnitBall_RK2
75 end subroutine
76
77 subroutine ndimInt2RK1()
79 call setLogVolUnitBall(logVolUnitBall_RK1, ndim_IK)
80 call getDummy_RK1()
81 end subroutine
82
83 subroutine ndimInt2RK2()
85 call setLogVolUnitBall(logVolUnitBall_RK2, ndim_IK)
86 call getDummy_RK2()
87 end subroutine
88
89 subroutine ndimRealRK1()
91 call setLogVolUnitBall(logVolUnitBall_RK1, ndim_RK1)
92 call getDummy_RK1()
93 end subroutine
94
95 subroutine ndimRealRK2()
97 call setLogVolUnitBall(logVolUnitBall_RK2, ndim_RK2)
98 call getDummy_RK2()
99 end subroutine
100
101end program benchmark
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Generate and return an object of type timing_type containing the benchmark timing information and sta...
Definition: pm_bench.F90:574
This module contains abstract interfaces and types that facilitate benchmarking of different procedur...
Definition: pm_bench.F90:41
integer, parameter RK2
Definition: pm_kind.F90:511
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 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 RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
integer, parameter RK1
Definition: pm_kind.F90:522
This is the class for creating benchmark and performance-profiling objects.
Definition: pm_bench.F90:386

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
7import os
8dirname = os.path.basename(os.getcwd())
9
10fontsize = 14
11
12df = pd.read_csv("main.out", delimiter = ",")
13colnames = list(df.columns.values)
14
15
18
19ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
20ax = plt.subplot()
21
22for colname in colnames[1:]:
23 plt.plot( df[colnames[0]].values
24 , df[colname].values
25 , linewidth = 2
26 )
27
28plt.xticks(fontsize = fontsize)
29plt.yticks(fontsize = fontsize)
30ax.set_xlabel(colnames[0], fontsize = fontsize)
31ax.set_ylabel("Runtime [ seconds ]", fontsize = fontsize)
32ax.set_title(" vs. ".join(colnames[1:])+"\nLower is better.", fontsize = fontsize)
33#ax.set_xscale("log")
34ax.set_yscale("log")
35plt.minorticks_on()
36plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
37ax.tick_params(axis = "y", which = "minor")
38ax.tick_params(axis = "x", which = "minor")
39ax.legend ( colnames[1:]
40 #, loc='center left'
41 #, bbox_to_anchor=(1, 0.5)
42 , fontsize = fontsize
43 )
44
45plt.tight_layout()
46plt.savefig("benchmark." + dirname + ".runtime.png")
47
48
51
52ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
53ax = plt.subplot()
54
55plt.plot( df[colnames[0]].values
56 , np.ones(len(df[colnames[0]].values))
57 , linestyle = "--"
58 #, color = "black"
59 , linewidth = 2
60 )
61for colname in colnames[2:]:
62 plt.plot( df[colnames[0]].values
63 , df[colname].values / df[colnames[1]].values
64 , linewidth = 2
65 )
66
67plt.xticks(fontsize = fontsize)
68plt.yticks(fontsize = fontsize)
69ax.set_xlabel(colnames[0], fontsize = fontsize)
70ax.set_ylabel("Runtime compared to {}".format(colnames[1]), fontsize = fontsize)
71ax.set_title("Runtime Ratio Comparison. Lower means faster.\nLower than 1 means faster than {}().".format(colnames[1]), fontsize = fontsize)
72#ax.set_xscale("log")
73#ax.set_yscale("log")
74plt.minorticks_on()
75plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
76ax.tick_params(axis = "y", which = "minor")
77ax.tick_params(axis = "x", which = "minor")
78ax.legend ( colnames[1:]
79 #, bbox_to_anchor = (1, 0.5)
80 #, loc = "center left"
81 , fontsize = fontsize
82 )
83
84plt.tight_layout()
85plt.savefig("benchmark." + dirname + ".runtime.ratio.png")

Visualization of the benchmark output

Benchmark moral
  1. The benchmark procedures named ndimInt2RK1 and ndimInt2RK2 call the generic interface setLogVolUnitBall with a ndim argument of type integer of kind IK and return the result as real of kind RKS and RKD respectively.
    Conversely, the procedures named ndimReal2RK1 and ndimReal2RK2 take ndim as real of kind RKS and RKD respectively and return results of the same type and kind as ndim.
    The first class of procedure interfaces (working with a ndim of type integer) compute the result using an iterative factorial computation approach.
    The second class of procedure interfaces (working with a ndim of type real) use the Fortran intrinsic log_gamma() to compute the results.
    This performance difference tends to be about a factor of 10 times for scalar optional arguments and grows substantially to larger factors with switching to increasing-size array-like optional dummy arguments.
  2. From the benchmark results it appears that the interface accepting ndim of type RK32 offers the fastest algorithm, followed by the iterative methods.
Test:
test_pm_ellipsoid


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, April 23, 2017, 1:36 AM, Institute for Computational Engineering and Sciences (ICES), University of Texas at Austin

Definition at line 558 of file pm_ellipsoid.F90.


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