ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation. |
This is the derived type for generating objects containing the specifications of a realization of an MMUVE distribution and a collection of randomly sampled points from it.
More...
Public Member Functions | |
procedure, pass | write => mmvue_type_write |
Public Attributes | |
type(range_type) | range |
The scalar object of type range_type containing the specifications of the MMVUE distribution. More... | |
integer(IK) | ndim |
The number of dimensions of the domain of the MMVUE distribution. More... | |
integer(IK) | nell |
The number of ellipsoids in the domain of the MMVUE distribution. More... | |
integer(IK) | nsam |
The number of random points sampled from the domain of the MMVUE distribution. More... | |
integer(IK) | nsim = 10000 |
The scalar of type integer of default kind IK, representing the number of simulated random vectors from the ellipsoids of the distribution to estimate the effective volumes of the ellipsoids (i.e., the inverse of the PDF of the distribution).More... | |
real(RKG) | logVolUnitBall |
The natural logarithm of the volume of a unit ball in ndim dimensions.More... | |
real(RKG) | logDensity |
The scalar of type real containing the natural logarithm of the inverse of the effective volume of a single sampled point from the distribution domain, normalized by the volume of a unit ball.In other words, logDensity = log(mmvue_type%nsam) - mmvue_type%logSumVolNormedEff .More... | |
real(RKG) | logSumVolNormedEff |
The natural logarithm of the sum of the effective volumes of the distribution ellipsoids, normalized by logVolUnitBall .This measure includes shared volumes between ellipsoids only once from one individual ellipsoid. Because this measure does not have closed expression, it must be approximated using Monte Carlo approach. By definition, this sum is always smaller than or equal to the sum of the individual normalized ellipsoid volumes. The equality happens only when no two ellipsoids overlap. More... | |
real(RKG), dimension(:,:), allocatable | std |
The matrix of shape (1:ndim, 1:nell) containing the standard deviations of the ellipsoids in the distribution.More... | |
real(RKG), dimension(:,:), allocatable | sample |
The matrix of shape (1:ndim, 1:nsam) containing the set of uniformly-sampled random points from the distribution domain.More... | |
real(RKG), dimension(:,:), allocatable | mean |
The matrix of shape (1:ndim, 1:nell) containing the set of nell ellipsoid centers of the distribution domain.More... | |
real(RKG), dimension(:,:,:), allocatable | invGram |
The matrix of shape (1:ndim, 1:ndim, 1:nell) containing the set of nell ellipsoid inverse-Gramian matrices of the distribution domain.More... | |
real(RKG), dimension(:,:,:), allocatable | choLowGramUpp |
The matrix of shape (1:ndim, 1:ndim + 1, 1:nell) containing the set of nell ellipsoid Cholesky-Gramian half (triangular) matrices of the distribution domain.More... | |
real(RKG), dimension(:), allocatable | cumPropVolNormed |
The vector of shape (0:nell) containing the cumulative proportion of volumes of individual ellipsoids in the distribution, used for generating the random sample .More... | |
real(RKG), dimension(:), allocatable | logVolNormed |
The vector of shape (1:nell) containing the natural logarithm of the set of nell volumes of the Gramian matrices of the distribution domain, normalized by the volume of a unit ball.More... | |
real(RKG), dimension(:,:), allocatable | mahalSq |
The matrix of shape (1:nell, 1:nsam) containing the square of the Mahalanobis distances the randomly sampled points from each ellipsoid center in the distribution.More... | |
real(RKG), dimension(:), allocatable | invmul |
The vector of shape (1:nsam) containing the inverse of sample multiplicity, that is, the ellipsoidal membership count of each randomly-sampled point from the distribution.For example, a multiplicity of 3 means that the isam sampled point falls within exactly three ellipsoids in the distribution.More... | |
integer(IK), dimension(:), allocatable | size |
The vector of shape (1:nell) containing the set of nell counts of randomly-sampled points from within each ellipsoid in the distribution.More... | |
integer(IK), dimension(:), allocatable | membership |
The vector of shape (1:nsam) containing the ellipsoidal membership of each randomly-sampled point from the distribution.More... | |
type(err_type) | err |
This is the derived type for generating objects containing the specifications of a realization of an MMUVE distribution and a collection of randomly sampled points from it.
See also the documentation of the type constructor mmvue_typer.
Possible calling interfaces ⛓
Example usage ⛓
ifort
compiler ⛓ ifort
compiler ⛓ gfortran
compiler ⛓ real
components of arbitrary kind, once GNU compiler support for PDTs is strong.
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.
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.
Definition at line 149 of file pm_clusTest.F90.
procedure, pass pm_clusTest::mmvue_type::write |
Definition at line 179 of file pm_clusTest.F90.
real(RKG), dimension(:,:,:), allocatable pm_clusTest::mmvue_type::choLowGramUpp |
The matrix of shape (1:ndim, 1:ndim + 1, 1:nell)
containing the set of nell
ellipsoid Cholesky-Gramian half (triangular) matrices of the distribution domain.
Definition at line 169 of file pm_clusTest.F90.
real(RKG), dimension(:), allocatable pm_clusTest::mmvue_type::cumPropVolNormed |
The vector of shape (0:nell)
containing the cumulative proportion of volumes of individual ellipsoids in the distribution, used for generating the random sample
.
Definition at line 170 of file pm_clusTest.F90.
type(err_type) pm_clusTest::mmvue_type::err |
Definition at line 177 of file pm_clusTest.F90.
real(RKG), dimension(:,:,:), allocatable pm_clusTest::mmvue_type::invGram |
The matrix of shape (1:ndim, 1:ndim, 1:nell)
containing the set of nell
ellipsoid inverse-Gramian matrices of the distribution domain.
Definition at line 168 of file pm_clusTest.F90.
real(RKG), dimension(:), allocatable pm_clusTest::mmvue_type::invmul |
The vector of shape (1:nsam)
containing the inverse of sample multiplicity, that is, the ellipsoidal membership count of each randomly-sampled point from the distribution.
For example, a multiplicity of 3
means that the isam
sampled point falls within exactly three ellipsoids in the distribution.
Definition at line 173 of file pm_clusTest.F90.
real(RKG) pm_clusTest::mmvue_type::logDensity |
The scalar of type real
containing the natural logarithm of the inverse of the effective volume of a single sampled point from the distribution domain, normalized by the volume of a unit ball.
In other words, logDensity = log(mmvue_type%nsam) - mmvue_type%logSumVolNormedEff
.
Definition at line 157 of file pm_clusTest.F90.
real(RKG) pm_clusTest::mmvue_type::logSumVolNormedEff |
The natural logarithm of the sum of the effective volumes of the distribution ellipsoids, normalized by logVolUnitBall
.
This measure includes shared volumes between ellipsoids only once from one individual ellipsoid.
Because this measure does not have closed expression, it must be approximated using Monte Carlo approach.
By definition, this sum is always smaller than or equal to the sum of the individual normalized ellipsoid volumes.
The equality happens only when no two ellipsoids overlap.
Definition at line 160 of file pm_clusTest.F90.
real(RKG), dimension(:), allocatable pm_clusTest::mmvue_type::logVolNormed |
The vector of shape (1:nell)
containing the natural logarithm of the set of nell
volumes of the Gramian matrices of the distribution domain, normalized by the volume of a unit ball.
Definition at line 171 of file pm_clusTest.F90.
real(RKG) pm_clusTest::mmvue_type::logVolUnitBall |
The natural logarithm of the volume of a unit ball in ndim
dimensions.
Definition at line 156 of file pm_clusTest.F90.
real(RKG), dimension(:,:), allocatable pm_clusTest::mmvue_type::mahalSq |
The matrix of shape (1:nell, 1:nsam)
containing the square of the Mahalanobis distances the randomly sampled points from each ellipsoid center in the distribution.
Definition at line 172 of file pm_clusTest.F90.
real(RKG), dimension(:,:), allocatable pm_clusTest::mmvue_type::mean |
The matrix of shape (1:ndim, 1:nell)
containing the set of nell
ellipsoid centers of the distribution domain.
Definition at line 167 of file pm_clusTest.F90.
integer(IK), dimension(:), allocatable pm_clusTest::mmvue_type::membership |
The vector of shape (1:nsam)
containing the ellipsoidal membership of each randomly-sampled point from the distribution.
Definition at line 176 of file pm_clusTest.F90.
integer(IK) pm_clusTest::mmvue_type::ndim |
The number of dimensions of the domain of the MMVUE distribution.
Definition at line 151 of file pm_clusTest.F90.
integer(IK) pm_clusTest::mmvue_type::nell |
The number of ellipsoids in the domain of the MMVUE distribution.
Definition at line 152 of file pm_clusTest.F90.
integer(IK) pm_clusTest::mmvue_type::nsam |
The number of random points sampled from the domain of the MMVUE distribution.
Definition at line 153 of file pm_clusTest.F90.
integer(IK) pm_clusTest::mmvue_type::nsim = 10000 |
The scalar of type integer
of default kind IK, representing the number of simulated random vectors from the ellipsoids of the distribution to estimate the effective volumes of the ellipsoids (i.e., the inverse of the PDF of the distribution).
Definition at line 154 of file pm_clusTest.F90.
type(range_type) pm_clusTest::mmvue_type::range |
The scalar object of type range_type containing the specifications of the MMVUE distribution.
Definition at line 150 of file pm_clusTest.F90.
real(RKG), dimension(:,:), allocatable pm_clusTest::mmvue_type::sample |
The matrix of shape (1:ndim, 1:nsam)
containing the set of uniformly-sampled random points from the distribution domain.
Definition at line 166 of file pm_clusTest.F90.
integer(IK), dimension(:), allocatable pm_clusTest::mmvue_type::size |
The vector of shape (1:nell)
containing the set of nell
counts of randomly-sampled points from within each ellipsoid in the distribution.
Definition at line 175 of file pm_clusTest.F90.
real(RKG), dimension(:,:), allocatable pm_clusTest::mmvue_type::std |
The matrix of shape (1:ndim, 1:nell)
containing the standard deviations of the ellipsoids in the distribution.
Definition at line 165 of file pm_clusTest.F90.