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

Generate and return the logarithm of the sum of the exponential of two input logarithmic values robustly (without causing overflow). More...

Detailed Description

Generate and return the logarithm of the sum of the exponential of two input logarithmic values robustly (without causing overflow).

Parameters
[in]smaller: The input scalar, or array of the same rank, shape, and size as other array-like input arguments, of either
  • type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128), or
  • type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
representing the smaller of the two logarithmic values to be added.
If it is complex, then the real component of the complex value must be smaller than the real component of the (larger) complex value.
(optional, it must be present if and only if the input argument larger is also present.)
[in]larger: The input scalar, or array of the same rank, shape, and size as other array-like input arguments, of the same type and kind as the input argument smaller, representing the larger of the two logarithmic values to be added.
If it is complex, then the real component of the complex value must be larger than the real component of the (larger) complex value.
(optional, it must be present if and only if the input argument smaller is also present.)
[in]minMax: The input vector of length 2 of either
  • type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128) or,
  • type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
representing the smaller and the larger of the two logarithmic values to be added as minMax(1) and minMax(2).
If it is complex, then the real component of the complex value must be smaller than the real component of the (larger) complex value.
If the smaller and larger values are not known a priori, the vector minMax can be instead constructed by calling getMinMax.
(optional, it must be present if and only if the input arguments smaller and larger are missing.)
control: The input scalar object that can be,
  1. the constant selection or equivalently, an object of type selection_type.
    Specifying this option enables the runtime checks for underflow occurrence via branching and dynamic dispatch.
    Enabling this option can aid runtime efficiency when the smaller to larger ratio is expected to cause underflow.
    In such cases, the exponentiation is avoided if control = selection, leading to faster runtime by avoiding exponentiation since it is highly expensive (on the order of ~200 CPU cycles).
    See the relevant benchmark here.
The presence of this argument is merely for compile-time resolution of the procedures of this generic interface.
(optional. If missing, a sequence control flow will be assumed without dynamic dispatch.)
Returns
logAddExp : The output scalar of the same type and kind as the input larger and smaller representing the logarithm of the sum of the exponential of two input logarithmic values computed robustly.


Possible calling interfaces

use pm_mathLogAddExp, only: getLogAddExp, selection
logAddExp = getLogAddExp(smaller, larger)
logAddExp = getLogAddExp(smaller, larger, control)
logAddExp = getLogAddExp(minMax(1:2))
logAddExp = getLogAddExp(minMax(1:2), control)
Generate and return the logarithm of the sum of the exponential of two input logarithmic values robus...
This module contains procedures and generic interfaces for adding two real or complex values without ...
Warning
As the argument names suggest, the value of the input argument larger must be larger than or equal to the value of the input argument smaller.
When the input arguments are of type complex, this condition must hold for the real components of the numbers.
These conditions are 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 this generic interface are elemental when the input arguments smaller and larger are present.
Note
The smaller and larger input values can be quickly obtained by a one-line call to setMinMax.
Alternatively, the smaller and larger values can be sorted in a minMax(1:2) array by calling getMinMax and passing the result directly to getLogAddExp.
See also
getLog1p
get1mexp
getMinMax
setMinMax
getLogAddExp
getLogSubExp
getLogSumExp


Example usage

1program example
2
3 use pm_kind, only: SK
4 use pm_io, only: display_type
5 use pm_kind, only: IK, CKS, CKD, CKH, RKS, RKD, RKH
7 use pm_distUnif, only: getUnifRand
8 use pm_mathMinMax, only: getMinMax
9
10 implicit none
11
12 type(display_type) :: disp
13
14 disp = display_type(file = "main.out.F90")
15
16 call disp%skip()
17 call disp%show("!%%%%%%%%%%%")
18 call disp%show("!32-bit real")
19 call disp%show("!%%%%%%%%%%%")
20 call disp%skip()
21
22 call disp%show("exp( getLogAddExp(smaller = log(4._RKS), larger = log(10._RKS)) )")
23 call disp%show( exp( getLogAddExp(smaller = log(4._RKS), larger = log(10._RKS)) ) )
24
25 call disp%skip()
26 call disp%show("!%%%%%%%%%%%")
27 call disp%show("!64-bit real")
28 call disp%show("!%%%%%%%%%%%")
29 call disp%skip()
30
31 call disp%show("exp( getLogAddExp(smaller = log(4._RKD), larger = log(10._RKD)) )")
32 call disp%show( exp( getLogAddExp(smaller = log(4._RKD), larger = log(10._RKD)) ) )
33
34 call disp%skip()
35 call disp%show("!%%%%%%%%%%%%")
36 call disp%show("!128-bit real")
37 call disp%show("!%%%%%%%%%%%%")
38 call disp%skip()
39
40 call disp%show("exp( getLogAddExp(smaller = log(4._RKH), larger = log(10._RKH)) )")
41 call disp%show( exp( getLogAddExp(smaller = log(4._RKH), larger = log(10._RKH)) ) )
42
43 call disp%skip()
44 call disp%show("!%%%%%%%%%%%%%%")
45 call disp%show("!32-bit complex")
46 call disp%show("!%%%%%%%%%%%%%%")
47 call disp%skip()
48
49 call disp%show("exp( getLogAddExp(smaller = log((4._CKS, -4._CKS)), larger = log((10._CKS, -10._CKS))) )")
50 call disp%show( exp( getLogAddExp(smaller = log((4._CKS, -4._CKS)), larger = log((10._CKS, -10._CKS))) ) )
51
52 call disp%skip()
53 call disp%show("!%%%%%%%%%%%%%%")
54 call disp%show("!64-bit complex")
55 call disp%show("!%%%%%%%%%%%%%%")
56 call disp%skip()
57
58 call disp%show("exp( getLogAddExp(smaller = log((4._CKD, -4._CKD)), larger = log((10._CKD, -10._CKD))) )")
59 call disp%show( exp( getLogAddExp(smaller = log((4._CKD, -4._CKD)), larger = log((10._CKD, -10._CKD))) ) )
60
61 call disp%skip()
62 call disp%show("!%%%%%%%%%%%%%%%")
63 call disp%show("!128-bit complex")
64 call disp%show("!%%%%%%%%%%%%%%%")
65 call disp%skip()
66
67 call disp%show("exp( getLogAddExp(smaller = log((4._CKH, -4._CKH)), larger = log((10._CKH, -10._CKH))) )")
68 call disp%show( exp( getLogAddExp(smaller = log((4._CKH, -4._CKH)), larger = log((10._CKH, -10._CKH))) ) )
69
70 call disp%skip()
71 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
72 call disp%show("!Input arguments as a vector of [minimum, maximum] by calling getMinMax()")
73 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
74 call disp%skip()
75
76 call disp%show("exp(getLogAddExp(getMinMax(getUnifRand(-1., 1., 2_IK))))")
77 call disp%show( exp(getLogAddExp(getMinMax(getUnifRand(-1., 1., 2_IK)))) )
78
79end program example
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
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 an array of size two containing the minimum and maximum of the two input values i...
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 CKH
The scalar integer constant of intrinsic default kind, representing the highest-precision complex kin...
Definition: pm_kind.F90:843
integer, parameter CKS
The single-precision complex kind in Fortran mode. On most platforms, this is a 32-bit real kind.
Definition: pm_kind.F90:570
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 CKD
The double precision complex kind in Fortran mode. On most platforms, this is a 64-bit real kind.
Definition: pm_kind.F90:571
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
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
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
This module contains procedures and generic interfaces for finding the minimum and maximum of two inp...
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!32-bit real
4!%%%%%%%%%%%
5
6exp( getLogAddExp(smaller = log(4._RKS), larger = log(10._RKS)) )
7+14.0000010
8
9!%%%%%%%%%%%
10!64-bit real
11!%%%%%%%%%%%
12
13exp( getLogAddExp(smaller = log(4._RKD), larger = log(10._RKD)) )
14+14.000000000000004
15
16!%%%%%%%%%%%%
17!128-bit real
18!%%%%%%%%%%%%
19
20exp( getLogAddExp(smaller = log(4._RKH), larger = log(10._RKH)) )
21+13.9999999999999999999999999999999985
22
23!%%%%%%%%%%%%%%
24!32-bit complex
25!%%%%%%%%%%%%%%
26
27exp( getLogAddExp(smaller = log((4._CKS, -4._CKS)), larger = log((10._CKS, -10._CKS))) )
28(+14.0000010, -14.0000010)
29
30!%%%%%%%%%%%%%%
31!64-bit complex
32!%%%%%%%%%%%%%%
33
34exp( getLogAddExp(smaller = log((4._CKD, -4._CKD)), larger = log((10._CKD, -10._CKD))) )
35(+13.999999999999998, -13.999999999999996)
36
37!%%%%%%%%%%%%%%%
38!128-bit complex
39!%%%%%%%%%%%%%%%
40
41exp( getLogAddExp(smaller = log((4._CKH, -4._CKH)), larger = log((10._CKH, -10._CKH))) )
42(+13.9999999999999999999999999999999985, -13.9999999999999999999999999999999969)
43
44!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
45!Input arguments as a vector of [minimum, maximum] by calling getMinMax()
46!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
47
48exp(getLogAddExp(getMinMax(getUnifRand(-1., 1., 2_IK))))
49+3.03160238
50
Benchmarks:


Benchmark :: The effects of control argument on runtime efficiency

The following program compares the runtime performance of getLogAddExp algorithm with and without checking for underflows.
1! Test the performance of `getLogAddExp()` with and without the `control` argument.
2program benchmark
3
4 use iso_fortran_env, only: error_unit
5 use pm_bench, only: bench_type
6 use pm_kind, only: IK, LK, RK, SK
7
8 implicit none
9
10 integer(IK) :: i
11 integer(IK) :: iarr
12 integer(IK) :: fileUnitN
13 integer(IK) :: fileUnitU
14 integer(IK) , parameter :: NARR = 11_IK
15 integer(IK) , parameter :: NBENCH = 2_IK
16 integer(IK) :: arraySize(NARR)
17 real(RK) :: dummySum = 0._RK
18 real(RK) :: maxArray
19 real(RK) , allocatable :: Larger(:)
20 real(RK) , allocatable :: Smaller(:)
21 real(RK) , allocatable :: LogAddExp(:)
22 type(bench_type) :: bench(2,NBENCH)
23 logical(LK) :: underflowEnabled
24
25 bench(1:2,1) = bench_type(name = SK_"getLogAddExpSequence", exec = getLogAddExpSequence , overhead = setOverhead)
26 bench(1:2,2) = bench_type(name = SK_"getLogAddExpSelection", exec = getLogAddExpSelection , overhead = setOverhead)
27
28 arraySize = [( 2_IK**iarr, iarr = 1_IK, NARR )]
29
30 write(*,"(*(g0,:,' '))")
31 write(*,"(*(g0,:,' '))") "getLogAddExp(...) vs. getLogAddExp(..., control) vs. direct method."
32 write(*,"(*(g0,:,' '))")
33
34 open(newunit = fileUnitN, file = "main.normal.out", status = "replace")
35 open(newunit = fileUnitU, file = "main.underflow.out", status = "replace")
36
37 write(fileUnitN, "(*(g0,:,','))") "arraySize", (bench(1,i)%name, i = 1, NBENCH)
38 write(fileUnitU, "(*(g0,:,','))") "arraySize", (bench(2,i)%name, i = 1, NBENCH)
39
40 loopOverArraySize: do iarr = 1, NARR
41
42 allocate(Larger(arraySize(iarr)))
43 allocate(Smaller(arraySize(iarr)))
44 allocate(LogAddExp(arraySize(iarr)))
45 write(*,"(*(g0,:,' '))") "Benchmarking with array size", arraySize(iarr)
46
47 underflowEnabled = .false._LK
48 do i = 1, NBENCH
49 bench(1,i)%timing = bench(1,i)%getTiming(minsec = 0.07_RK)
50 end do
51 write(fileUnitN,"(*(g0,:,','))") arraySize(iarr), (bench(1,i)%timing%mean, i = 1, NBENCH)
52
53 underflowEnabled = .true._LK
54 do i = 1, NBENCH
55 bench(2,i)%timing = bench(2,i)%getTiming(minsec = 0.07_RK)
56 end do
57 write(fileUnitU,"(*(g0,:,','))") arraySize(iarr), (bench(2,i)%timing%mean, i = 1, NBENCH)
58
59 deallocate(LogAddExp, Smaller, Larger)
60
61 end do loopOverArraySize
62 write(*,"(*(g0,:,' '))") dummySum
63 write(*,"(*(g0,:,' '))")
64
65 close(fileUnitN)
66 close(fileUnitU)
67
68contains
69
70 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
71 ! procedure wrappers.
72 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
73
74 subroutine setOverhead()
75 call setArray()
76 call getDummy()
77 end subroutine
78
79 subroutine setArray()
80 use pm_mathMinMax, only: setMinMax
81 call random_number(Larger)
82 call random_number(Smaller)
83 if (underflowEnabled) then
84 Larger = Larger * (maxexponent(0._RK) - minexponent(0._RK)) + minexponent(0._RK)
85 Smaller = Smaller * (maxexponent(0._RK) - minexponent(0._RK)) + minexponent(0._RK)
86 end if
87 call setMinMax(Smaller, Larger)
88 end subroutine
89
90 subroutine getDummy()
91 dummySum = dummySum + LogAddExp(1)
92 end subroutine
93
94 subroutine getLogAddExpSequence()
96 call setArray()
97 LogAddExp(:) = getLogAddExp(Smaller, Larger)
98 call getDummy()
99 end subroutine
100
101 subroutine getLogAddExpSelection()
102 use pm_mathLogAddExp, only: getLogAddExp, selection
103 call setArray()
104 LogAddExp(:) = getLogAddExp(Smaller, Larger, selection)
105 call getDummy()
106 end subroutine
107
108end program benchmark
Generate and return an object of type timing_type containing the benchmark timing information and sta...
Definition: pm_bench.F90:574
Return the minimum and maximum of the two input scalar values a and b in a and b, respectively,...
This module contains abstract interfaces and types that facilitate benchmarking of different procedur...
Definition: pm_bench.F90:41
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 LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
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
12
15
16df = pd.read_csv("main.normal.out")
17colnames = list(df.columns.values)
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 #, np.ones(len(df[colnames[0]].values))
25 , df[colname].values / df[colnames[1]].values
26 , linewidth = 2
27 #, linestyle = "-"
28 )
29
30plt.xticks(fontsize = fontsize)
31plt.yticks(fontsize = fontsize)
32ax.set_xlabel("Array Size", fontsize = fontsize)
33ax.set_ylabel("Runtime Ratio ( W.R.T. cenabled = .false. )", fontsize = fontsize)
34ax.set_title("Performance when the input arrays\ncause no underflow instances.\nLower is better.", fontsize = fontsize)
35ax.set_xscale("log")
36#ax.set_yscale("log")
37plt.minorticks_on()
38plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
39ax.tick_params(axis = "y", which = "minor")
40ax.tick_params(axis = "x", which = "minor")
41ax.legend ( colnames[1:]
42 #, loc='center left'
43 #, bbox_to_anchor=(1, 0.5)
44 , fontsize = fontsize
45 )
46
47plt.tight_layout()
48plt.savefig("benchmark." + dirname + ".normal.png")
49
50
53
54df = pd.read_csv("main.underflow.out")
55colnames = list(df.columns.values)
56
57ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
58ax = plt.subplot()
59
60for colname in colnames[1:]:
61 plt.plot( df[colnames[0]].values
62 #, np.ones(len(df[colnames[0]].values))
63 , df[colname].values / df[colnames[1]].values
64 , linewidth = 2
65 #, linestyle = "-"
66 )
67
68plt.xticks(fontsize = fontsize)
69plt.yticks(fontsize = fontsize)
70ax.set_xlabel("Array Size", fontsize = fontsize)
71ax.set_ylabel("Runtime Ratio W.R.T. {})".format(colnames[1]), fontsize = fontsize)
72ax.set_title("Performance when the input arrays\ncause many underflow instances.\nLower is better.", fontsize = fontsize)
73ax.set_xscale("log")
74#ax.set_yscale("log")
75plt.minorticks_on()
76plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
77ax.tick_params(axis = "y", which = "minor")
78ax.tick_params(axis = "x", which = "minor")
79ax.legend ( colnames[1:]
80 #, loc='center left'
81 #, bbox_to_anchor=(1, 0.5)
82 , fontsize = fontsize
83 )
84
85plt.tight_layout()
86plt.savefig("benchmark." + dirname + ".underflow.png")

Visualization of the benchmark output

Benchmark moral
  1. If the ratio of the smaller to larger input arguments causes frequent underflows, then setting control = selection to allow branching will likely result in a faster runtime.
    Conversely, if the divisions are not expected to cause any or too many underflows, then removing the input argument control when calling getLogSubExp will lead to a better runtime performance (at the expense of occasional expensive but redundant exponentiation operations).
  2. If the procedure getLogAddExp is to be called billions of times, then it would make sense to manually inline the procedure implementation in your code as procedure call will have a non-negligible performance overhead.
Test:
test_pm_mathLogAddExp


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, Thursday 1:45 AM, August 22, 2019, Dallas, TX

Definition at line 166 of file pm_mathLogAddExp.F90.


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