s_nlp_hkmeans_core Subroutine

public pure subroutine s_nlp_hkmeans_core(x, nd, nv, nc, gm, cm, cl, cc, cov, sigma)

Perform agglomerative hierarchical clustering using centroid linkage and the Mahalanobis distance, then passes cluster centroids and covariance matrix to kmeans cluster procedure for refinement.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x(nd,nv)

input data matrix (samples, variables)

integer(kind=i4), intent(in) :: nd

number of data points

integer(kind=i4), intent(in) :: nv

number of variables

integer(kind=i4), intent(in) :: nc

number of clusters (target)

real(kind=wp), intent(out) :: gm(nv)

global means for each variable

real(kind=wp), intent(out) :: cm(nv,nc)

cluster centroids

integer(kind=i4), intent(out) :: cl(nd)

cluster assignments for each data point

integer(kind=i4), intent(out) :: cc(nc)

cluster sizes

real(kind=wp), intent(out) :: cov(nv,nv)

covariance matrix

real(kind=wp), intent(out) :: sigma(nv)

standard deviation per variable


Calls

proc~~s_nlp_hkmeans_core~~CallsGraph proc~s_nlp_hkmeans_core s_nlp_hkmeans_core proc~s_nlp_hclust_core s_nlp_hclust_core proc~s_nlp_hkmeans_core->proc~s_nlp_hclust_core proc~s_nlp_kmeans_core s_nlp_kmeans_core proc~s_nlp_hkmeans_core->proc~s_nlp_kmeans_core proc~f_lin_mahalanobis_core f_lin_mahalanobis_core proc~s_nlp_hclust_core->proc~f_lin_mahalanobis_core proc~f_sts_cov_core f_sts_cov_core proc~s_nlp_hclust_core->proc~f_sts_cov_core proc~f_sts_mean_core f_sts_mean_core proc~s_nlp_hclust_core->proc~f_sts_mean_core proc~f_sts_var_core f_sts_var_core proc~s_nlp_hclust_core->proc~f_sts_var_core proc~s_utl_sort s_utl_sort proc~s_nlp_hclust_core->proc~s_utl_sort proc~s_nlp_kmeans_core->proc~f_lin_mahalanobis_core proc~s_nlp_kmeans_core->proc~f_sts_cov_core proc~s_nlp_kmeans_core->proc~f_sts_mean_core proc~s_nlp_kmeans_core->proc~f_sts_var_core proc~s_nlp_kmeans_core->proc~s_utl_sort proc~f_lin_mahalanobis_core->proc~f_sts_cov_core proc~s_utl_cholesky_solve s_utl_cholesky_solve proc~f_lin_mahalanobis_core->proc~s_utl_cholesky_solve proc~f_sts_cov_core->proc~f_sts_mean_core proc~f_sts_var_core->proc~f_sts_mean_core chol chol proc~s_utl_cholesky_solve->chol

Called by

proc~~s_nlp_hkmeans_core~~CalledByGraph proc~s_nlp_hkmeans_core s_nlp_hkmeans_core proc~s_nlp_hkmeans s_nlp_hkmeans proc~s_nlp_hkmeans->proc~s_nlp_hkmeans_core interface~fsml_hkmeans fsml_hkmeans interface~fsml_hkmeans->proc~s_nlp_hkmeans

Source Code

pure subroutine s_nlp_hkmeans_core(x, nd, nv, nc, gm, cm, cl, cc, cov, sigma)

! ==== Description
!! Perform agglomerative hierarchical clustering using centroid linkage
!! and the Mahalanobis distance, then passes cluster centroids and
!! covariance matrix to kmeans cluster procedure for refinement.

! ==== Declarations
  real(wp)   , intent(in)  :: x(nd, nv)     !! input data matrix (samples, variables)
  integer(i4), intent(in)  :: nd            !! number of data points
  integer(i4), intent(in)  :: nv            !! number of variables
  integer(i4), intent(in)  :: nc            !! number of clusters (target)
  real(wp)   , intent(out) :: gm(nv)        !! global means for each variable
  real(wp)   , intent(out) :: cm(nv,nc)     !! cluster centroids
  integer(i4), intent(out) :: cl(nd)        !! cluster assignments for each data point
  integer(i4), intent(out) :: cc(nc)        !! cluster sizes
  real(wp)   , intent(out) :: cov(nv,nv)    !! covariance matrix
  real(wp)   , intent(out) :: sigma(nv)     !! standard deviation per variable
  real(wp)                 :: cm_h(nv,nc)   !! cluster centroids passed from hclust to kmeans
  real(wp)                 :: cov_h(nv,nv)  !! covariance matrix passed from hclust to kmeans

! ==== Instructions

  ! hierarchical clustering
  call s_nlp_hclust_core(x, nd, nv, nc, gm, cm_h, cl, cc, cov_h, sigma)

  ! kmeans refinement
  call s_nlp_kmeans_core(x, nd, nv, nc, cm_h, gm, cm, cl, cc, &
                       & cov, sigma, cov_in=cov_h)

end subroutine s_nlp_hkmeans_core