ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_clusTest::mmvue_type Type Reference

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_typerange
 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_typeerr
 

Detailed Description

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

type(mmvue_type) :: mmvue
mmvue = mmvue_type(rng, range = range_type(), nsim = nsim)
This module contains procedures and routines for the creating test datasets for clustering algorithms...
Definition: pm_clusTest.F90:39
This is the derived type for generating objects containing the specifications of a realization of an ...
This is the derived type for generating objects containing the range of specifications of an MMVUE di...
Definition: pm_clusTest.F90:79
See also
mmvue_typer


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_clusTest, only: mmvue_type, rngf_type
5 use pm_io, only: getErrTableWrite, trans
6 use pm_io, only: display_type
7 use pm_val2str, only: getStr
8
9 implicit none
10
11 integer(IK) :: itry
12 type(mmvue_type) :: mmvue
13 type(display_type) :: disp
14 disp = display_type(file = "main.out.F90")
15
16 do itry = 0, 9
17 call disp%skip()
18 call disp%show("mmvue = mmvue_type(rng = rngf_type())")
19 mmvue = mmvue_type(rng = rngf_type())
20 call disp%show("[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]")
21 call disp%show( [mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim] )
22 call disp%show("if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'")
23 if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
24 call disp%skip()
25 end do
26
27end program example
Generate and return the iostat code resulting from writing the input table of rank 1 or 2 to the spec...
Definition: pm_io.F90:5940
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 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
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
2mmvue = mmvue_type(rng = rngf_type())
3[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
4+2, +11, +5000, +10000
5if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
6
7
8mmvue = mmvue_type(rng = rngf_type())
9[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
10+2, +13, +5000, +10000
11if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
12
13
14mmvue = mmvue_type(rng = rngf_type())
15[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
16+2, +7, +5000, +10000
17if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
18
19
20mmvue = mmvue_type(rng = rngf_type())
21[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
22+2, +4, +5000, +10000
23if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
24
25
26mmvue = mmvue_type(rng = rngf_type())
27[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
28+2, +5, +5000, +10000
29if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
30
31
32mmvue = mmvue_type(rng = rngf_type())
33[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
34+2, +17, +5000, +10000
35if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
36
37
38mmvue = mmvue_type(rng = rngf_type())
39[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
40+2, +13, +5000, +10000
41if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
42
43
44mmvue = mmvue_type(rng = rngf_type())
45[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
46+2, +17, +5000, +10000
47if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
48
49
50mmvue = mmvue_type(rng = rngf_type())
51[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
52+2, +1, +5000, +10000
53if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
54
55
56mmvue = mmvue_type(rng = rngf_type())
57[mmvue%ndim, mmvue%nell, mmvue%nsam, mmvue%nsim]
58+2, +1, +5000, +10000
59if (0 /= getErrTableWrite(SK_'mmvue_type.'//getStr(itry)//SK_'.txt', reshape([transpose(mmvue%sample), real(mmvue%membership, kind(mmvue%sample))], [mmvue%nsam, mmvue%ndim + 1_IK]))) error stop 'Failed to write the random vectors file.'
60
61

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
8import os
9procname = os.path.basename(os.getcwd())
10
11linewidth = 2
12fontsize = 17
13
14pattern = procname + "*.txt"
15fileList = glob.glob(pattern)
16
17for 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, 8), dpi = 200)
28
29
30 plt.rcParams.update({'font.size': fontsize - 2})
31 ax = plt.subplot()
32 #ax = fig.add_axes([left, bottom, width, height]) # scatter plot
33 #ax_histx = fig.add_axes([left, bottom + height + spacing, width, 0.2], sharex = ax) # histx
34 #ax_histy = fig.add_axes([left + width + spacing, bottom, 0.2, height], sharey = ax) # histy
35
36 for axes in [ax]:#, ax_histx, ax_histy]:
37 axes.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
38 axes.tick_params(axis = "y", which = "minor")
39 axes.tick_params(axis = "x", which = "minor")
40
41 # no labels
42 #ax_histy.tick_params(axis = "y", labelleft = False)
43 #ax_histx.tick_params(axis = "x", labelbottom = False)
44
45 # the scatter plot:
46 ax.scatter ( df.values[:, 0]
47 , df.values[:, 1]
48 , c = df.values[:, -1]
49 , cmap = "gist_rainbow"
50 , zorder = 1000
51 , marker = "o"
52 , alpha = 1
53 , s = 0.7
54 )
55
56 #ax_histx.hist(df.values[:, 0], bins = 50, zorder = 1000)
57 #ax_histy.hist(df.values[:, 1], bins = 50, orientation = "horizontal", zorder = 1000)
58
59 ax.set_xlabel("X", fontsize = 17)
60 ax.set_ylabel("Y", fontsize = 17)
61 #plt.title("file: " + file)
62
63 plt.savefig(file.replace(".txt",".svg"))

Visualization of the example output



















Test:
test_pm_clusTest
Todo:
High Priority: This type must be generalized to accept 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.

  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, 2012, 12:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 149 of file pm_clusTest.F90.

Member Function/Subroutine Documentation

◆ write()

procedure, pass pm_clusTest::mmvue_type::write

Definition at line 179 of file pm_clusTest.F90.

Member Data Documentation

◆ choLowGramUpp

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.

◆ cumPropVolNormed

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.

◆ err

type(err_type) pm_clusTest::mmvue_type::err

Definition at line 177 of file pm_clusTest.F90.

◆ invGram

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.

◆ invmul

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.

◆ logDensity

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.

◆ logSumVolNormedEff

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.

◆ logVolNormed

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.

◆ logVolUnitBall

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.

◆ mahalSq

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.

◆ mean

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.

◆ membership

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.

◆ ndim

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.

◆ nell

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.

◆ nsam

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.

◆ nsim

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.

◆ range

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.

◆ sample

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.

◆ size

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.

◆ std

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.


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