Generate and return the cosmological size (or equivalently, the age) of the Universe at the desired redshift normalized to Hubble Distance (or equivalently, to Hubble Time), 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 Universe Size at a given redshift \(z\) is simply related to the Hubble Parameter as,
\begin{equation}
\large
D_S(z; \Omega_M, \Omega_\Lambda, \Omega_R, \Omega_K) =
D_H \int_{z}^{+\infty} \frac{1}{(1+z') E(z'; \Omega_M, \Omega_\Lambda, \Omega_R, \Omega_K)} ~ dz' ~,
\end{equation}
or equivalently, the Universe Age to a cosmological object as redshift \(z\) is simply related to the Hubble Parameter as,
\begin{equation}
\large
T_A(z; \Omega_M, \Omega_\Lambda, \Omega_R, \Omega_K) =
T_H \int_{z}^{+\infty} \frac{1}{(1+z') E(z'; \Omega_M, \Omega_\Lambda, \Omega_R, \Omega_K)} ~ dz' ~,
\end{equation}
where,
- \(z\) is the redshift,
- \(T_A\) is the cosmological Universe Age,
- \(D_S\) is the cosmological Universe Size,
- \(E(z; \cdots)\) is the dimensionless Hubble Parameter,
- \(D_H = \frac{C}{H_0}\) is the Hubble Distance,
- \(T_H = \frac{1}{H_0}\) is the Hubble Time,
- \(H_0\) is the Hubble Constant.
- \(C\) is the speed of light,
The Universe Size, also known as the Light-Travel Size, is defined as the distance traveled by light since the Big Bang to a given redshift.
Dividing the Universe Size by the Speed of Light yields the Universe Age.
For instance, the radius of the observable universe in this distance measure becomes the age of the universe multiplied by the speed of light (1 light year/year), which turns out to be approximately 13.8
billion light years.
The value returned by the procedures under this generic interface is \(\frac{D_S}{D_H}\).
The value returned by the procedures under this generic interface is \(\frac{T_A}{T_H}\).
The default method of computing the Universe Size is numerical integration via Adaptive Global Gauss-Kronrod 10-21 Quadrature rule.
- 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 Universe Size 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] | 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
sizeUnivNormed
: 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 Universe Size (or Universe Age) at the desired redshift normalized to the Hubble Distance (or Hubble Time).
Possible calling interfaces ⛓
sizeUnivNormed
= getSizeUnivNormed(zplus1, omegaM, omegaL, reltol, neval
= neval,
err = err)
sizeUnivNormed
= getSizeUnivNormed(zplus1, omegaM, omegaL, omegaR, reltol, neval
= neval,
err = err)
sizeUnivNormed
= getSizeUnivNormed(zplus1, omegaM, omegaL, omegaR, omegaK, sqrtAbsOmegaK, reltol, neval
= neval,
err = err)
Generate and return the cosmological size (or equivalently, the age) of the Universe at the desired r...
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
.
-
Note that for some cosmological model parameters, the integral in the definition of the Universe Size/Age is divergent.
- 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
getDisLookbackNormed
getSizeUnivNormed
getDisAngNormed
getDisLumNormed
getDisComNormed
LOG_HUBBLE_CONST
HUBBLE_DISTANCE_MPC
HUBBLE_CONST
OMEGA_M
OMEGA_L
OMEGA_R
OMEGA_K
Example usage ⛓
10 type(display_type) :: disp
15 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
16 call disp%show(
"!Compute the Comoving diameter distance in units of Hubble Distance.")
17 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
21 call disp%show(
"getSizeUnivNormed(zplus1 = 1., reltol = sqrt(epsilon(0.)))")
26 call disp%show(
"getSizeUnivNormed(zplus1 = 1.1, reltol = sqrt(epsilon(0.)))")
31 call disp%show(
"getSizeUnivNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.6, reltol = 0.0001)")
36 call disp%show(
"getSizeUnivNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.4, omegaR = 0.2, reltol = 0.0001)")
41 call disp%show(
"getSizeUnivNormed(zplus1 = [ 2., 3.], omegaM = 0.4, omegaL = 0.4, omegaR = 0.4, omegaK = -0.2, reltol = 0.0001)")
42 call disp%show(
getSizeUnivNormed(zplus1
= [
2.,
3.], omegaM
= 0.4, omegaL
= 0.4, omegaR
= 0.4, omegaK
= -0.2, reltol
= 0.0001) )
53 real,
allocatable :: zplus1(:), OmegaM(:), OmegaL(:)
54 integer :: fileUnit, i
56 OmegaM
= [
0.2,
0.3,
1.0]
58 zplus1
= 1. + getLogSpace(
log(
0.0001),
log(
10000.),
500_IK)
60 open(newunit
= fileUnit, file
= "getSizeUnivNormed.RK.txt")
61 write(fileUnit,
"(*(g0,:,','))")
"z", (
"DisLum_"//getStr(OmegaM(i),
"(g0.1)")
//"_"//getStr(OmegaL(i),
"(g0.1)"), i
= 1,
size(OmegaM))
62 do i
= 1,
size(zplus1)
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.
This is a generic method of the derived type display_type with pass attribute.
Generate and return the conversion of the input value to an output Fortran string,...
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
real(RKB), parameter HUBBLE_DISTANCE_GLY
The scalar real constant of kind with highest available precision RKB representing the Hubble Distanc...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
This module contains the generic procedures for converting values of different types and kinds to For...
Generate and return an object of type display_type.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
16+0.361917049,
+0.201027289
19getSizeUnivNormed(zplus1
= [
2.,
3.], omegaM
= 0.4, omegaL
= 0.4, omegaR
= 0.2, reltol
= 0.0001)
20+0.216050461,
+0.103596434
23getSizeUnivNormed(zplus1
= [
2.,
3.], omegaM
= 0.4, omegaL
= 0.4, omegaR
= 0.4, omegaK
= -0.2, reltol
= 0.0001)
24+0.174306065,
+0.803307965E-1
Postprocessing of the example output ⛓
3import matplotlib.pyplot
as plt
15xlab = {
"CK" :
"redshift: Z ( real/imaginary components )"
16 ,
"IK" :
"redshift: Z ( integer-valued )"
17 ,
"RK" :
"redshift: Z ( real-valued )"
19legends = [
"$\Omega_M = 0.2, \Omega_\Lambda = 0.8$"
20 ,
"$\Omega_M = 0.3, \Omega_\Lambda = 0.7$"
21 ,
"$\Omega_M = 1.0, \Omega_\Lambda = 0.0$"
24for kind
in [
"IK",
"CK",
"RK"]:
26 pattern =
"*." + kind +
".txt"
27 fileList = glob.glob(pattern)
28 if len(fileList) == 1:
30 df = pd.read_csv(fileList[0], delimiter =
",")
32 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
36 plt.plot( df.values[:, 0]
37 , df.values[:,1:len(legends)+1]
41 plt.plot( df.values[:, 1]
42 , df.values[:,1:len(legends)+1]
47 plt.plot( df.values[:, 0]
48 , df.values[:,1:len(legends)+1]
56 plt.xticks(fontsize = fontsize - 2)
57 plt.yticks(fontsize = fontsize - 2)
58 ax.set_xlabel(xlab[kind], fontsize = 17)
59 ax.set_ylabel(
"Universe Size [ GLY ]", 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")
67 plt.savefig(fileList[0].replace(
".txt",
".png"))
69 elif len(fileList) > 1:
71 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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- Author:
- Amir Shahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
Definition at line 270 of file pm_cosmology.F90.