SF_STAT
Description
SciFortran module for statistical estimators and visualization
Quick access
- Types:
- Routines:
pdf_allocate()
,pdf_deallocate()
,pdf_save()
,pdf_read()
,pdf_set_range()
,pdf_push_sigma()
,pdf_sigma()
,pdf_accumulate()
,pdf_normalize()
,pdf_print()
,pdf_mean()
,pdf_var()
,pdf_sdev()
,pdf_moment()
,pdf_skew()
,pdf_curt()
,pdf_print_moments()
,get_moments()
,get_mean()
,get_sd()
,get_var()
,get_skew()
,get_curt()
,get_covariance()
,histogram_allocate()
,histogram_deallocate()
,histogram_reset()
,histogram_set_range_uniform()
,histogram_accumulate()
,histogram_get_range()
,histogram_get_value()
,histogram_print()
Used modules
sf_arrays
: SciFortran module for array creation and manipulationsf_integrate
: SciFortran module for function integrationsf_iotools
: SciFortran module for reading and writing datasf_linalg
: SciFortran module for linear algebra
Types
- type sf_stat/histogram
- Type fields:
n [integer]
range (•) [real, pointer]
bin (•) [real, pointer]
- type sf_stat/pdf_kernel
- Type fields:
n [integer]
xmin [real]
xmax [real]
dx [real]
ndata [integer]
x (•) [real, allocatable]
pdf (•) [real, allocatable]
sigma [real]
status [logical]
variance [logical]
rescale [logical]
- type sf_stat/pdf_kernel_2d
- Type fields:
n (2) [integer]
xmin (2) [real]
xmax (2) [real]
dx (2) [real]
ndata [integer]
x (•) [real, allocatable]
y (•) [real, allocatable]
pdf (•, •) [real, allocatable]
sigma (2, 2) [real]
status [logical]
variance [logical]
rescale [logical]
Subroutines and functions
- interface sf_stat/pdf_allocate(self, n, nvec)
- Parameters:
self [pdf_kernel]
n [integer]
nvec (2) [integer]
- interface sf_stat/pdf_deallocate(self)
- Parameters:
self [pdf_kernel]
- interface sf_stat/pdf_save(self, pfile)
- Parameters:
self [pdf_kernel]
pfile [character(len=*)]
- interface sf_stat/pdf_read(self, pfile)
- Parameters:
self [pdf_kernel]
pfile [character(len=*)]
- interface sf_stat/pdf_set_range(self, a, b)
- Parameters:
self [pdf_kernel]
a (various shapes) [real]
b (various shapes) [real]
- interface sf_stat/pdf_push_sigma(self, sigma)
- Parameters:
self [pdf_kernel]
sigma (various shapes) [real]
- interface sf_stat/pdf_sigma(self, data, h, sdev, n, nvec)
- Parameters:
self [pdf_kernel]
data (various shapes) [real]
h (various shapes) [real] – Silverman’s rule of thumb.
sdev (various shapes) [real] – Silverman’s rule of thumb.
n [integer]
nvec (2) [integer]
- interface sf_stat/pdf_accumulate(self, data[, sigma])
- Parameters:
self [pdf_kernel]
data (various shapes) [real]
- Options:
sigma (various shapes) [real]
- interface sf_stat/pdf_normalize(self)
- Parameters:
self [pdf_kernel]
- interface sf_stat/pdf_print(self, pfile[, normalize])
- Parameters:
self [pdf_kernel]
pfile [character(len=*)]
- Options:
normalize [logical]
- interface sf_stat/pdf_mean(self)
- Parameters:
self [pdf_kernel]
- Return:
mean [real]
- interface sf_stat/pdf_var(self)
- Parameters:
self [pdf_kernel]
- Return:
var [real]
- interface sf_stat/pdf_sdev(self)
- Parameters:
self [pdf_kernel]
- Return:
sdev [real]
- interface sf_stat/pdf_moment(self, n[, mu])
- Parameters:
self [pdf_kernel]
n [integer]
- Options:
mu [real]
- Return:
mom [real]
- interface sf_stat/pdf_skew(self)
- Parameters:
self [pdf_kernel]
- Return:
skew [real]
- interface sf_stat/pdf_curt(self)
- Parameters:
self [pdf_kernel]
- Return:
curt [real]
- interface sf_stat/pdf_print_moments(self, pfile)
- Parameters:
self [pdf_kernel]
pfile [character(len=*)]
- subroutine sf_stat/get_moments(data[, ave, sdev, var, skew, curt])
- Parameters:
data (•) [real, in]
- Options:
ave [real, out]
sdev [real, out]
var [real, out]
skew [real, out]
curt [real, out]
- function sf_stat/get_mean(data)
- Parameters:
data (•) [real, in]
- Return:
mean [real]
- function sf_stat/get_sd(data)
- Parameters:
data (•) [real, in]
- Return:
sd [real]
- function sf_stat/get_var(data)
- Parameters:
data (•) [real, in]
- Return:
var [real]
- function sf_stat/get_skew(data)
- Parameters:
data (•) [real, in]
- Return:
skew [real]
- function sf_stat/get_curt(data)
- Parameters:
data (•) [real, in]
- Return:
curt [real]
- function sf_stat/get_covariance(data, mean)
- Parameters:
data (•, •) [real, in,required]
mean (size(data, 1)) [real, in]
- Return:
covariance (size(data, 1), size(data, 1)) [real]
- subroutine sf_stat/histogram_set_range_uniform(h, xmin, xmax)
- Parameters:
h [histogram, inout]
xmin [real, in]
xmax [real, in]
- subroutine sf_stat/histogram_accumulate(h, x, w)
- Parameters:
h [histogram, inout]
x [real, in]
w [real, in]
- subroutine sf_stat/histogram_get_range(h, index, lower, upper)
- Parameters:
h [histogram, in]
index [integer, in]
lower [real, out]
upper [real, out]