X-Git-Url: https://git.auder.net/?a=blobdiff_plain;f=epclust%2FR%2Fclustering.R;h=3993e7685c97b194644d3600ceeea6b7bdac54ac;hb=24ed5d835e2eebaaa4d5f8296f8d2e2132cc6398;hp=6090517c6b6464d4c253ba52b8efdf29cb56c823;hpb=0e2dce80a3fddaca50c96c6c27a8b32468095d6c;p=epclust.git diff --git a/epclust/R/clustering.R b/epclust/R/clustering.R index 6090517..3993e76 100644 --- a/epclust/R/clustering.R +++ b/epclust/R/clustering.R @@ -1,70 +1,196 @@ -# Cluster one full task (nb_curves / ntasks series) -clusteringTask = function(indices, ncores) +#' @name clustering +#' @rdname clustering +#' @aliases clusteringTask1 computeClusters1 computeClusters2 +#' +#' @title Two-stage clustering, withing one task (see \code{claws()}) +#' +#' @description \code{clusteringTask1()} runs one full stage-1 task, which consists in +#' iterated stage 1 clustering (on nb_curves / ntasks energy contributions, computed +#' through discrete wavelets coefficients). \code{computeClusters1()} and +#' \code{computeClusters2()} correspond to the atomic clustering procedures respectively +#' for stage 1 and 2. The former applies the clustering algorithm (PAM) on a +#' contributions matrix, while the latter clusters a chunk of series inside one task +#' (~max nb_series_per_chunk) +#' +#' @param indices Range of series indices to cluster in parallel (initial data) +#' @param getContribs Function to retrieve contributions from initial series indices: +#' \code{getContribs(indices)} outpus a contributions matrix +#' @param contribs matrix of contributions (e.g. output of \code{curvesToContribs()}) +#' @inheritParams computeSynchrones +#' @inheritParams claws +#' +#' @return For \code{clusteringTask1()} and \code{computeClusters1()}, the indices of the +#' computed (K1) medoids. Indices are irrelevant for stage 2 clustering, thus +#' \code{computeClusters2()} outputs a big.matrix of medoids +#' (of size limited by nb_series_per_chunk) +NULL + +#' @rdname clustering +#' @export +clusteringTask1 = function( + indices, getContribs, K1, nb_series_per_chunk, ncores_clust=1, verbose=FALSE, parll=TRUE) { - cl = parallel::makeCluster(ncores) - parallel::clusterExport(cl, - varlist=c("K1","getCoefs"), - envir=environment()) - repeat + if (verbose) + cat(paste("*** Clustering task on ",length(indices)," lines\n", sep="")) + + wrapComputeClusters1 = function(inds) { + if (parll) + require("epclust", quietly=TRUE) + if (verbose) + cat(paste(" computeClusters1() on ",length(inds)," lines\n", sep="")) + inds[ computeClusters1(getContribs(inds), K1) ] + } + + if (parll) { - nb_workers = max( 1, round( length(indices_clust) / nb_series_per_chunk ) ) - indices_workers = lapply(seq_len(nb_workers), function(i) { - upper_bound = ifelse( i K1) + { + indices_workers = .spreadIndices(indices, nb_series_per_chunk) + if (parll) + indices = unlist( parallel::parLapply(cl, indices_workers, wrapComputeClusters1) ) + else + indices = unlist( lapply(indices_workers, wrapComputeClusters1) ) + } + if (parll) + parallel::stopCluster(cl) + + indices #medoids } -# Apply the clustering algorithm (PAM) on a coeffs or distances matrix -computeClusters1 = function(indices, getCoefs, K1) - indices[ cluster::pam(getCoefs(indices), K1, diss=FALSE)$id.med ] +#' @rdname clustering +#' @export +computeClusters1 = function(contribs, K1) + cluster::pam(contribs, K1, diss=FALSE)$id.med -# Cluster a chunk of series inside one task (~max nb_series_per_chunk) -computeClusters2 = function(indices, K2, getSeries, to_file) +#' @rdname clustering +#' @export +computeClusters2 = function(medoids, K2, + getRefSeries, nb_ref_curves, nb_series_per_chunk, ncores_clust=1,verbose=FALSE,parll=TRUE) { - if (is.null(indices)) + synchrones = computeSynchrones(medoids, + getRefSeries, nb_ref_curves, nb_series_per_chunk, ncores_clust, verbose, parll) + distances = computeWerDists(synchrones, ncores_clust, verbose, parll) + #TODO: if PAM cannot take big.matrix in input, cast it before... (more than OK in RAM) + medoids[ cluster::pam(distances, K2, diss=TRUE)$medoids , ] +} + +#' computeSynchrones +#' +#' Compute the synchrones curves (sum of clusters elements) from a matrix of medoids, +#' using L2 distances. +#' +#' @param medoids big.matrix of medoids (curves of same length as initial series) +#' @param getRefSeries Function to retrieve initial series (e.g. in stage 2 after series +#' have been replaced by stage-1 medoids) +#' @param nb_ref_curves How many reference series? (This number is known at this stage) +#' @inheritParams claws +#' +#' @return A big.matrix of size K1 x L where L = data_length +#' +#' @export +computeSynchrones = function(medoids, getRefSeries, + nb_ref_curves, nb_series_per_chunk, ncores_clust=1,verbose=FALSE,parll=TRUE) +{ + + + +#TODO: si parll, getMedoids + serialization, pass only getMedoids to nodes +# --> BOF... chaque node chargera tous les medoids (efficacité) :/ ==> faut que ça tienne en RAM +#au pire :: C-ifier et charger medoids 1 by 1... + + #MIEUX :: medoids DOIT etre une big.matrix partagée ! + + computeSynchronesChunk = function(indices) { - #get series from file + if (verbose) + cat(paste("--- Compute synchrones for ",length(indices)," lines\n", sep="")) + ref_series = getRefSeries(indices) + #get medoids indices for this chunk of series + for (i in seq_len(nrow(ref_series))) + { + j = which.min( rowSums( sweep(medoids, 2, ref_series[i,], '-')^2 ) ) + if (parll) + synchronicity::lock(m) + synchrones[j,] = synchrones[j,] + ref_series[i,] + counts[j,1] = counts[j,1] + 1 + if (parll) + synchronicity::unlock(m) + } } -#Puis K-means après WER... - if (WER=="mix" > 0) + + K = nrow(medoids) + # Use bigmemory (shared==TRUE by default) + synchronicity to fill synchrones in // + # TODO: if size > RAM (not our case), use file-backed big.matrix + synchrones = bigmemory::big.matrix(nrow=K,ncol=ncol(medoids),type="double",init=0.) + counts = bigmemory::big.matrix(nrow=K,ncol=1,type="double",init=0) + # synchronicity is only for Linux & MacOS; on Windows: run sequentially + parll = (requireNamespace("synchronicity",quietly=TRUE) + && parll && Sys.info()['sysname'] != "Windows") + if (parll) + m <- synchronicity::boost.mutex() + + if (parll) { - curves = computeSynchrones(indices) - dists = computeWerDists(curves) - indices = computeClusters(dists, K2, diss=TRUE) + cl = parallel::makeCluster(ncores_clust) + parallel::clusterExport(cl, + varlist=c("synchrones","counts","verbose","medoids","getRefSeries"), + envir=environment()) } - if (to_file) - #write results to file (JUST series ; no possible ID here) -} -# Compute the synchrones curves (sum of clusters elements) from a clustering result -computeSynchrones = function(inds) - sapply(seq_along(inds), colMeans(getSeries(inds[[i]]$indices,inds[[i]]$ids))) + indices_workers = .spreadIndices(seq_len(nb_ref_curves), nb_series_per_chunk) + ignored <- + if (parll) + parallel::parLapply(indices_workers, computeSynchronesChunk) + else + lapply(indices_workers, computeSynchronesChunk) + + if (parll) + parallel::stopCluster(cl) + + #TODO: can we avoid this loop? ( synchrones = sweep(synchrones, 1, counts, '/') ) + for (i in seq_len(K)) + synchrones[i,] = synchrones[i,] / counts[i,1] + #NOTE: odds for some clusters to be empty? (when series already come from stage 2) + # ...maybe; but let's hope resulting K1' be still quite bigger than K2 + noNA_rows = sapply(seq_len(K), function(i) all(!is.nan(synchrones[i,]))) + if (all(noNA_rows)) + return (synchrones) + # Else: some clusters are empty, need to slice synchrones + synchrones[noNA_rows,] +} -# Compute the WER distance between the synchrones curves (in columns) -computeWerDist = function(curves) +#' computeWerDists +#' +#' Compute the WER distances between the synchrones curves (in rows), which are +#' returned (e.g.) by \code{computeSynchrones()} +#' +#' @param synchrones A big.matrix of synchrones, in rows. The series have same length +#' as the series in the initial dataset +#' @inheritParams claws +#' +#' @return A big.matrix of size K1 x K1 +#' +#' @export +computeWerDists = function(synchrones, ncores_clust=1,verbose=FALSE,parll=TRUE) { - if (!require("Rwave", quietly=TRUE)) - stop("Unable to load Rwave library") - n <- nrow(curves) - delta <- ncol(curves) + + + +#TODO: re-organize to call computeWerDist(x,y) [C] (in //?) from two indices + big.matrix + + + n <- nrow(synchrones) + delta <- ncol(synchrones) #TODO: automatic tune of all these parameters ? (for other users) nvoice <- 4 - # noctave = 2^13 = 8192 half hours ~ 180 days ; ~log2(ncol(curves)) + # noctave = 2^13 = 8192 half hours ~ 180 days ; ~log2(ncol(synchrones)) noctave = 13 # 4 here represent 2^5 = 32 half-hours ~ 1 day #NOTE: default scalevector == 2^(0:(noctave * nvoice) / nvoice) * s0 (?) - scalevector <- 2^(4:(noctave * nvoice) / nvoice) * 2 + scalevector <- 2^(4:(noctave * nvoice) / nvoice + 1) #condition: ( log2(s0*w0/(2*pi)) - 1 ) * nvoice + 1.5 >= 1 s0=2 w0=2*pi @@ -72,33 +198,106 @@ computeWerDist = function(curves) s0log = as.integer( (log2( s0*w0/(2*pi) ) - 1) * nvoice + 1.5 ) totnoct = noctave + as.integer(s0log/nvoice) + 1 - # (normalized) observations node with CWT - Xcwt4 <- lapply(seq_len(n), function(i) { - ts <- scale(ts(curves[,i]), center=TRUE, scale=scaled) - totts.cwt = Rwave::cwt(ts,totnoct,nvoice,w0,plot=0) + computeCWT = function(i) + { + if (verbose) + cat(paste("+++ Compute Rwave::cwt() on serie ",i,"\n", sep="")) + ts <- scale(ts(synchrones[i,]), center=TRUE, scale=scaled) + totts.cwt = Rwave::cwt(ts, totnoct, nvoice, w0, plot=FALSE) ts.cwt = totts.cwt[,s0log:(s0log+noctave*nvoice)] #Normalization sqs <- sqrt(2^(0:(noctave*nvoice)/nvoice)*s0) - sqres <- sweep(ts.cwt,MARGIN=2,sqs,'*') + sqres <- sweep(ts.cwt,2,sqs,'*') sqres / max(Mod(sqres)) - }) + } - Xwer_dist <- matrix(0., n, n) + if (parll) + { + cl = parallel::makeCluster(ncores_clust) + parallel::clusterExport(cl, + varlist=c("synchrones","totnoct","nvoice","w0","s0log","noctave","s0","verbose"), + envir=environment()) + } + + # list of CWT from synchrones + # TODO: fit in RAM, OK? If not, 2 options: serialize, compute individual distances + Xcwt4 <- + if (parll) + parallel::parLapply(cl, seq_len(n), computeCWT) + else + lapply(seq_len(n), computeCWT) + + if (parll) + parallel::stopCluster(cl) + + Xwer_dist <- bigmemory::big.matrix(nrow=n, ncol=n, type="double") fcoefs = rep(1/3, 3) #moving average on 3 values (TODO: very slow! correct?!) - for (i in 1:(n-1)) + if (verbose) + cat("*** Compute WER distances from CWT\n") + + #TODO: computeDistances(i,j), et répartir les n(n-1)/2 couples d'indices + #là c'est trop déséquilibré + + computeDistancesLineI = function(i) { + if (verbose) + cat(paste(" Line ",i,"\n", sep="")) for (j in (i+1):n) { - #TODO: later, compute CWT here (because not enough storage space for 200k series) - # 'circular=TRUE' is wrong, should just take values on the sides; to rewrite in C + #TODO: 'circular=TRUE' is wrong, should just take values on the sides; to rewrite in C num <- filter(Mod(Xcwt4[[i]] * Conj(Xcwt4[[j]])), fcoefs, circular=TRUE) WX <- filter(Mod(Xcwt4[[i]] * Conj(Xcwt4[[i]])), fcoefs, circular=TRUE) WY <- filter(Mod(Xcwt4[[j]] * Conj(Xcwt4[[j]])), fcoefs, circular=TRUE) - wer2 <- sum(colSums(num)^2) / sum( sum(colSums(WX) * colSums(WY)) ) + wer2 <- sum(colSums(num)^2) / sum( sum(colSums(WX) * colSums(WY)) ) + if (parll) + synchronicity::lock(m) Xwer_dist[i,j] <- sqrt(delta * ncol(Xcwt4[[1]]) * (1 - wer2)) Xwer_dist[j,i] <- Xwer_dist[i,j] + if (parll) + synchronicity::unlock(m) } + Xwer_dist[i,i] = 0. } - diag(Xwer_dist) <- numeric(n) + + parll = (requireNamespace("synchronicity",quietly=TRUE) + && parll && Sys.info()['sysname'] != "Windows") + if (parll) + m <- synchronicity::boost.mutex() + + ignored <- + if (parll) + { + parallel::mclapply(seq_len(n-1), computeDistancesLineI, + mc.cores=ncores_clust, mc.allow.recursive=FALSE) + } + else + lapply(seq_len(n-1), computeDistancesLineI) + Xwer_dist[n,n] = 0. Xwer_dist } + +# Helper function to divide indices into balanced sets +.spreadIndices = function(indices, nb_per_chunk) +{ + L = length(indices) + nb_workers = floor( L / nb_per_chunk ) + if (nb_workers == 0) + { + # L < nb_series_per_chunk, simple case + indices_workers = list(indices) + } + else + { + indices_workers = lapply( seq_len(nb_workers), function(i) + indices[(nb_per_chunk*(i-1)+1):(nb_per_chunk*i)] ) + # Spread the remaining load among the workers + rem = L %% nb_per_chunk + while (rem > 0) + { + index = rem%%nb_workers + 1 + indices_workers[[index]] = c(indices_workers[[index]], indices[L-rem+1]) + rem = rem - 1 + } + } + indices_workers +}