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

Generate and return the cosmological Luminosity Distance normalized to Hubble Distance, given the user-specified cosmological parameters. More...

Detailed Description

Generate and return the cosmological Luminosity Distance normalized to Hubble Distance, given the user-specified cosmological parameters.

Assuming \((\Omega_M, \Omega_\Lambda, \Omega_R, \Omega_K)\) represent the normalized densities of Dark Matter, Dark Energy, Radiation Energy, and Curvature in a Universe with negligible neutrino mass such that \(\Omega_M + \Omega_\Lambda + \Omega_R + \Omega_K = 1\), the Luminosity Distance to a cosmological object as redshift \(z\) is simply related to the Transverse Comoving Distance as,

\begin{eqnarray} \large D_L(z; \Omega_M, \Omega_\Lambda, \Omega_R, \Omega_K) &=& (z + 1) ~ D_M(z; \Omega_M, \Omega_\Lambda, \Omega_R, \Omega_K) ~, \\ &=& (z + 1)^2 ~ D_A(z; \Omega_M, \Omega_\Lambda, \Omega_R, \Omega_K) ~, \end{eqnarray}

where,

  1. \(D_L(z; \cdots)\) is the cosmological Luminosity Distance in units of Mpc as a function of redshift,
  2. \(D_M(z; \cdots)\) is the cosmological Transverse Comoving Distance in units of Mpc as a function of redshift,
  3. \(D_A(z; \cdots)\) is the cosmological Angular Diameter Distance in units of Mpc as a function of redshift,
  4. \(z\) is the redshift,

The value returned by the procedures under this generic interface is \(\frac{D_L}{D_H}\).
The default method of computing the Luminosity Distance is the same as that of the Transverse Comoving Distance.

Parameters
[in]zplus1: The input scalar or array of the same rank as other array-like arguments, of type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128) representing the redshift plus one, \(\log(z+1)\), at which the Luminosity Distance must be computed.
[in]omegaM: The input scalar or array of the same rank as other array-like arguments, of the same type and kind as the input argument zplus1 representing the normalized matter density in the universe.
(optional, default = OMEGA_M. It must be present if omegaL ( \(\Omega_\Lambda\)) is also present.)
[in]omegaL: The input scalar or array of the same rank as other array-like arguments, of the same type and kind as the input argument zplus1 representing the normalized Dark Energy density in the universe.
(optional, default = OMEGA_L. It must be present if omegaR is also present.)
[in]omegaR: The input scalar or array of the same rank as other array-like arguments, of the same type and kind as the input argument zplus1 representing the normalized radiation density in the universe.
(optional, default = OMEGA_R. It must be present if omegaK is also present.)
[in]omegaK: The input scalar or array of the same rank as other array-like arguments, of the same type and kind as the input argument zplus1 representing the normalized curvature density of the universe.
(optional, default = OMEGA_K)
[in]sqrtAbsOmegaK: The input scalar or array of the same rank as other array-like arguments, of the same type and kind as the input argument zplus1 representing the square root of the normalized curvature density of the universe sqrt(omegaK).
This argument is requested as an input to avoid the redundant costly sqrt(OmegaK) operation within the algorithm when the procedure is to be called repeatedly with the same omegaK.
(optional, it must be present if and only if omegaK is also present.)
[in]reltol: See the description of the corresponding argument in the documentation of getQuadErr.
A reasonable recommended value is reltol = sqrt(epsilon(real(0, kind(zplus1)))).
[out]neval: See the description of the corresponding argument in the documentation of getQuadErr.
(optional. If missing, the number of function evaluations will not be returned.)
[out]err: The output scalar of type integer of default kind IK, that is set to zero if the integration converges without any errors.
Otherwise, a non-zero value of err indicates the occurrence of an error of varying severities.
See the description of the corresponding output argument in the documentation of getQuadErr.
Returns
disLumNormed : The output scalar or array of the same rank as other array-like arguments, of the same type and kind as the input argument zplus1 containing the cosmological Luminosity Distance at the desired redshift normalized to the Hubble Distance.


Possible calling interfaces

disLumNormed = getDisLumNormed(zplus1, reltol, neval = neval, err = err)
disLumNormed = getDisLumNormed(zplus1, omegaM, omegaL, reltol, neval = neval, err = err)
disLumNormed = getDisLumNormed(zplus1, omegaM, omegaL, omegaR, reltol, neval = neval, err = err)
disLumNormed = getDisLumNormed(zplus1, omegaM, omegaL, omegaR, omegaK, sqrtAbsOmegaK, reltol, neval = neval, err = err)
Generate and return the cosmological Luminosity Distance normalized to Hubble Distance,...
This module contains procedures and generic interfaces and constants for cosmological calculations.
Warning
The conditions listed in the documentation of getDisComTransNormed must hold for this generic interface.
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 under discussion are elemental.
Note
Dropping all optional arguments corresponds to the \(\Lambda\)CDM Universe with the latest experimental parameter inferences.
Setting omegaM = 1 and omegaL = 0 corresponds to the Einstein–de Sitter model of the universe proposed by Albert Einstein and Willem de Sitter in 1932.
See also
getHubbleParamNormedSq
getDisComTransNormed
getDisAngNormed
getDisLumNormed
getDisComNormed
LOG_HUBBLE_CONST
HUBBLE_DISTANCE_MPC
HUBBLE_CONST
OMEGA_M
OMEGA_L
OMEGA_R
OMEGA_K


Example usage

1program example
2
3 use pm_kind, only: SK, IK
5 use pm_io, only: display_type
6
7 implicit none
8
9 type(display_type) :: disp
10
11 disp = display_type(file = "main.out.F90")
12
13 call disp%skip()
14 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
15 call disp%show("!Compute the luminosity distance in units of Hubble Distance.")
16 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
17 call disp%skip()
18
19 call disp%skip()
20 call disp%show("getDisLumNormed(zplus1 = 1., reltol = sqrt(epsilon(0.)))")
21 call disp%show( getDisLumNormed(zplus1 = 1., reltol = sqrt(epsilon(0.))) )
22 call disp%skip()
23
24 call disp%skip()
25 call disp%show("getDisLumNormed(zplus1 = 1.1, reltol = sqrt(epsilon(0.)))")
26 call disp%show( getDisLumNormed(zplus1 = 1.1, reltol = sqrt(epsilon(0.))) )
27 call disp%skip()
28
29 call disp%skip()
30 call disp%show("getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.6, reltol = 0.0001)")
31 call disp%show( getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.6, reltol = 0.0001) )
32 call disp%skip()
33
34 call disp%skip()
35 call disp%show("getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.4, omegaR = 0.2, reltol = 0.0001)")
36 call disp%show( getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.4, omegaR = 0.2, reltol = 0.0001) )
37 call disp%skip()
38
39 call disp%skip()
40 call disp%show("getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.4, omegaR = 0.4, omegaK = -0.2, sqrtAbsOmegaK = sqrt(0.2), reltol = 0.0001)")
41 call disp%show( getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.4, omegaR = 0.4, omegaK = -0.2, sqrtAbsOmegaK = sqrt(0.2), reltol = 0.0001) )
42 call disp%skip()
43
44 ! Generate both the cosmic rate and the rate density.
45
46 block
47
48 use pm_val2str, only: getStr
49 use pm_arraySpace, only: getLinSpace
50 use pm_arraySpace, only: getLogSpace
51
52 real, allocatable :: zplus1(:), OmegaM(:), OmegaL(:)
53 integer :: fileUnit, i
54
55 OmegaM = [0., 0.5, 1.]
56 OmegaL = 1. - OmegaM
57 zplus1 = 1. + getLogSpace(log(0.0001), log(10000.), 500_IK)
58
59 open(newunit = fileUnit, file = "getDisLumNormed.RK.txt")
60 write(fileUnit, "(*(g0,:,','))") "z", ("DisLum_"//getStr(OmegaM(i),"(g0.1)")//"_"//getStr(OmegaL(i),"(g0.1)"), i = 1, size(OmegaM))
61 do i = 1, size(zplus1)
62 write(fileUnit, "(*(g0,:,','))") zplus1(i) - 1., getDisLumNormed(zplus1(i), OmegaM, OmegaL, sqrt(epsilon(0.)))
63 end do
64 close(fileUnit)
65
66 end block
67
68end program example
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Generate count evenly-logarithmically-spaced points over the interval [base**logx1,...
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 conversion of the input value to an output Fortran string,...
Definition: pm_val2str.F90:167
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 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 module contains the generic procedures for converting values of different types and kinds to For...
Definition: pm_val2str.F90:58
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!Compute the luminosity distance in units of Hubble Distance.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7getDisLumNormed(zplus1 = 1., reltol = sqrt(epsilon(0.)))
8+0.00000000
9
10
11getDisLumNormed(zplus1 = 1.1, reltol = sqrt(epsilon(0.)))
12+0.107397571
13
14
15getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.6, reltol = 0.0001)
16+1.46069443, +3.35959578
17
18
19getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.4, omegaR = 0.2, reltol = 0.0001)
20+1.26433861, +2.70926046
21
22
23getDisLumNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.4, omegaR = 0.4, omegaK = -0.2, sqrtAbsOmegaK = sqrt(0.2), reltol = 0.0001)
24+1.16154969, +2.38689947
25
26

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
9fontsize = 17
10
11marker ={ "CK" : "-"
12 , "IK" : "."
13 , "RK" : "-"
14 }
15xlab = { "CK" : "redshift: Z ( real/imaginary components )"
16 , "IK" : "redshift: Z ( integer-valued )"
17 , "RK" : "redshift: Z ( real-valued )"
18 }
19legends = [ "$\Omega_M = 0.0, \Omega_\Lambda = 1.0$"
20 , "$\Omega_M = 0.5, \Omega_\Lambda = 0.5$"
21 , "$\Omega_M = 1.0, \Omega_\Lambda = 0.0$"
22 ]
23
24for kind in ["IK", "CK", "RK"]:
25
26 pattern = "*." + kind + ".txt"
27 fileList = glob.glob(pattern)
28 if len(fileList) == 1:
29
30 df = pd.read_csv(fileList[0], delimiter = ",")
31
32 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
33 ax = plt.subplot()
34
35 if kind == "CK":
36 plt.plot( df.values[:, 0]
37 , df.values[:,1:len(legends)+1]
38 , marker[kind]
39 #, color = "r"
40 )
41 plt.plot( df.values[:, 1]
42 , df.values[:,1:len(legends)+1]
43 , marker[kind]
44 #, color = "blue"
45 )
46 else:
47 plt.plot( df.values[:, 0]
48 , df.values[:,1:len(legends)+1]
49 , marker[kind]
50 #, color = "r"
51 )
52 ax.legend ( legends
53 , fontsize = fontsize
54 )
55
56 plt.xticks(fontsize = fontsize - 2)
57 plt.yticks(fontsize = fontsize - 2)
58 ax.set_xlabel(xlab[kind], fontsize = 17)
59 ax.set_ylabel("Luminosity Distance / Hubble Distance", fontsize = 17)
60 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
61 ax.tick_params(axis = "y", which = "minor")
62 ax.tick_params(axis = "x", which = "minor")
63 ax.set_xscale("log")
64 ax.set_yscale("log")
65
66 plt.savefig(fileList[0].replace(".txt",".png"))
67
68 elif len(fileList) > 1:
69
70 sys.exit("Ambiguous file list exists.")

Visualization of the example output
Test:
test_pm_cosmology


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, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin

Definition at line 2344 of file pm_cosmology.F90.


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